main.cpp 14.9 KB
Newer Older
Kirill Terekhov's avatar
Kirill Terekhov committed
1
2
#include <string>
#include <iostream>
3
4
#include <iomanip>
#include <cmath>
Kirill Terekhov's avatar
Kirill Terekhov committed
5

Kirill Terekhov's avatar
Kirill Terekhov committed
6
#include "inmost.h"
7
#include "inner_parser.h"
Kirill Terekhov's avatar
Kirill Terekhov committed
8
using namespace INMOST;
Igor Konshin's avatar
Igor Konshin committed
9

Kirill Terekhov's avatar
Kirill Terekhov committed
10
11
12
13
14
15
#if defined(USE_MPI)
#define BARRIER MPI_Barrier(MPI_COMM_WORLD);
#else
#define BARRIER
#endif

16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
int main(int argc, char ** argv) {
    int processRank = 0, processorsCount = 1;

    #if defined(USE_MPI)
    MPI_Init(&argc, &argv);
    MPI_Comm_rank(MPI_COMM_WORLD, &processRank);  // Get the rank of the current process
    MPI_Comm_size(MPI_COMM_WORLD, &processorsCount ); // Get the total number of processors used
    #endif

    if (processRank == 0) {
        std::cout << "Starting MatSolve2" << std::endl;
    }

    {
        std::string matrixFileName = "";
        std::string vectorBFileName = "";
        std::string vectorXFileName = "";
        std::string parametersFileName = "";
        int solverType = 0;

        bool matrixFound = false;
        bool vectorBFound = false;
        bool vectorXFound = false;
        bool parametersFound = false;
        bool typeFound = false;

        //Parse argv parameters
        if (argc == 1) goto helpMessage;
        int i;
        for (i = 1; i < argc; i++) {
            //Print help message and exit
            if (strcmp(argv[i], "-h") == 0 || strcmp(argv[i], "--help") == 0) {
                helpMessage:
                if (processRank == 0) {
                    std::cout << "Help message: " << std::endl;
                    std::cout << "Command line options: " << std::endl;
                    std::cout << "Required: " << std::endl;
                    std::cout << "-m, --matrix <Matrix file name>" << std::endl;
                    std::cout << "Optional: " << std::endl;
                    std::cout << "-b, --bvector <RHS vector file name>" << std::endl;
                    std::cout << "-x, --xvector <X vector file name>" << std::endl;
                    std::cout << "-p, --parameters <Solver parameters file name>" << std::endl;
                    std::cout << "-t, --type <Solver type index>" << std::endl;
                    std::cout << "    0: INNER_ILU2 " << std::endl;
                    std::cout << "    1: INNER_DDPQILUC " << std::endl;
                    std::cout << "    2: INNER_MPTILUC " << std::endl;
                    std::cout << "    3: INNER_MPTILU2 " << std::endl;
                    std::cout << "    4: Trilinos_Aztec " << std::endl;
                    std::cout << "    5: Trilinos_Belos " << std::endl;
                    std::cout << "    6: Trilinos_ML " << std::endl;
                    std::cout << "    7: Trilinos_Ifpack " << std::endl;
                    std::cout << "    8: PETSc " << std::endl;
                    std::cout << "    9: ANI " << std::endl;
                    std::cout << "    10: FCBIILU2 " << std::endl;
                    std::cout << "    11: K3BIILU2 " << std::endl;
                    std::cout << "    12: SUPERLU " << std::endl << std::endl;
                    std::cout << "-h, --help - print this message." << std::endl;
                }
                #if defined(USE_MPI)
                MPI_Finalize();
                #endif
                return 0;
            }
            //Matrix file name found with -m or --matrix options
            if (strcmp(argv[i], "-m") == 0 || strcmp(argv[i], "--matrix") == 0) {
                if (processRank == 0) {
                    std::cout << "Matrix file found: " << argv[i + 1] << std::endl;
                }
                matrixFound = true;
                matrixFileName = std::string(argv[i + 1]);
                i++;
                continue;
            }
            //B vector file name found with -b or --bvector options
            if (strcmp(argv[i], "-b") == 0 || strcmp(argv[i], "--bvector") == 0) {
                if (processRank == 0) {
                    std::cout << "B vector file found: " << argv[i + 1] << std::endl;
                }
                vectorBFound = true;
                vectorBFileName = std::string(argv[i + 1]);
                i++;
                continue;
            }
            //X vector file name found with -x or --xvector options
            if (strcmp(argv[i], "-x") == 0 || strcmp(argv[i], "--xvector") == 0) {
                if (processRank == 0) {
                    std::cout << "X vector file found: " << argv[i + 1] << std::endl;
                }
                vectorXFound = true;
                vectorXFileName = std::string(argv[i + 1]);
                i++;
                continue;
            }
            //Parameters file name found with -p or --parameters options
            if (strcmp(argv[i], "-p") == 0 || strcmp(argv[i], "--parameters") == 0) {
                if (processRank == 0) {
                    std::cout << "Solver parameters file found: " << argv[i + 1] << std::endl;
                }
                parametersFound = true;
                parametersFileName = std::string(argv[i + 1]);
                i++;
                continue;
            }
            //Solver type found with -t ot --type options
            if (strcmp(argv[i], "-t") == 0 || strcmp(argv[i], "--type") == 0) {
                if (processRank == 0) {
                    std::cout << "Solver type index found: " << argv[i + 1] << std::endl;
                }
                typeFound = true;
                solverType = atoi(argv[i + 1]);
                i++;
                continue;
            }
        }

        if (!matrixFound) {
            if (processRank == 0) {
                std::cout <<
                "Matrix not found, you can specify matrix file name using -m or --matrix options, otherwise specify -h option to see all options, exiting...";
            }
            return -1;
        }

        if (!typeFound) {
            if (processRank == 0) {
                std::cout <<
                "Solver type not found in command line, you can specify solver type with -t or --type option, using INNER_ILU2 solver by default." <<
                std::endl;
            }
        }

        if (!vectorBFound) {
            if (processRank == 0) {
                std::cout <<
                "B vector not found, you can specify b vector file name with -b or --bvector option, using identity vector by default." <<
                std::endl;
            }
        }


        Solver::Type type;
        switch (solverType) {
            case 0: type = Solver::INNER_ILU2; break;
            case 1: type = Solver::INNER_DDPQILUC; break;
            case 2: type = Solver::INNER_MPTILUC; break;
            case 3: type = Solver::INNER_MPTILU2; break;
            case 4: type = Solver::Trilinos_Aztec; break;
            case 5: type = Solver::Trilinos_Belos; break;
            case 6: type = Solver::Trilinos_ML; break;
            case 7: type = Solver::Trilinos_Ifpack; break;
            case 8: type = Solver::PETSc; break;
            case 9: type = Solver::ANI; break;
            case 10: type = Solver::FCBIILU2; break;
            case 11: type = Solver::K3BIILU2; break;
            case 12: type = Solver::SUPERLU; break;
            default:
                if (processRank == 0) {
                    std::cout << "Invalid solver type index: " << solverType << " , using INNER_ILU2 instead." <<
                    std::endl;
                }
                type = Solver::INNER_ILU2;
                break;
        }

        Solver::Initialize(&argc, &argv, parametersFound ? parametersFileName.c_str()
                                                         : NULL); // Initialize the linear solver in accordance with args

        if (processRank == 0) {
            std::cout << "Solving with " << Solver::TypeName(type) << std::endl;
        }

        Sparse::Matrix mat("A"); // Declare the matrix of the linear system to be solved
        Sparse::Vector b("rhs"); // Declare the right-hand side vector
        Sparse::Vector x("sol"); // Declare the solution vector
        double tempTimer = Timer(), solvingTimer;
        mat.Load(matrixFileName); //if interval parameters not set, matrix will be divided automatically
        BARRIER
        if (processRank == 0) {
            std::cout << "load matrix:    " << Timer() - tempTimer << std::endl;
        }
        tempTimer = Timer();
        if (vectorBFound) {
            b.Load(vectorBFileName); // Load RHS vector
            //b.Save(vectorBFileName + "_saved_test.rhs");
        } else { // Set local RHS to 1 if it was not specified
            INMOST_DATA_ENUM_TYPE mbeg, mend, k;
            mat.GetInterval(mbeg, mend);
            b.SetInterval(mbeg, mend);
            for (k = mbeg; k < mend; ++k) b[k] = 1.0;
        }
        BARRIER
        if (processRank == 0) {
            std::cout << "load vector:    " << Timer() - tempTimer << std::endl;
        }
        solvingTimer = Timer();
        int iters;
        bool success;
        double resid, realresid = 0;
        std::string reason;
        Solver s(type); // Declare the linear solver by specified type

        if (parametersFound) {
            char *fileName = findInnerOptions(parametersFileName.c_str());
            if (fileName != NULL) {
                InnerOptions *options = parseInnerDatabaseOptions(fileName);
221
222
                if (options != NULL) {
                    for (unsigned ii = 0; ii < options->options.size(); ii++) {
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
                        InnerOption *option = options->options[ii];
                        if (option->type == ENUM) {
                            s.SetParameterEnum(option->name, (unsigned int) atoi(option->value.c_str()));
                        } else {
                            s.SetParameterReal(option->name, atof(option->value.c_str()));
                        };
                    }
                    delete options;
                }
                free(fileName);
            }
        }
        BARRIER

        //s.SetParameterEnum("maximum_iterations", 1000);
        //s.SetParameterEnum("gmres_substeps", 4);
        //s.SetParameterReal("relative_tolerance", 1.0e-6);
        //s.SetParameterReal("absolute_tolerance", 1.0e-16);
        //s.SetParameterReal("divergence_tolerance", 1e+200);

        //s.SetParameterEnum("reorder_nonzeros", 0);
        //s.SetParameterEnum("rescale_iterations", 8);
        //s.SetParameterEnum("adapt_ddpq_tolerance", 0);

        //s.SetParameterReal("drop_tolerance", 3.0e-3);
        //s.SetParameterReal("reuse_tolerance", 1.0e-5);
        //s.SetParameterReal("ddpq_tolerance", 0.0);

        //s.SetParameterEnum("condition_estimation", 1);
        //s.SetParameterEnum("schwartz_overlap", 3);


        tempTimer = Timer();
        s.SetMatrix(mat); // Compute the preconditioner for the original matrix
        BARRIER
        if (processRank == 0) std::cout << "preconditioner time: " << Timer() - tempTimer << std::endl;
        tempTimer = Timer();
        success = s.Solve(b, x); // Solve the linear system with the previously computted preconditioner
        BARRIER
        solvingTimer = Timer() - solvingTimer;
        if (processRank == 0) std::cout << "iterations time:     " << Timer() - tempTimer << std::endl;
        iters = s.Iterations(); // Get the number of iterations performed
        resid = s.Residual();   // Get the final residual achieved
        reason = s.GetReason(); // Get the convergence reason
        //x.Save("output.sol");  // Save the solution if required

        // Compute the true residual
        double aresid = 0, bresid = 0;
        Sparse::Vector test;
        tempTimer = Timer();
        Solver::OrderInfo info;
        info.PrepareMatrix(mat, 0);
        info.PrepareVector(x);
        info.Update(x);

        mat.MatVec(1.0, x, 0.0, test); // Multiply the original matrix by a vector
        {
            INMOST_DATA_ENUM_TYPE mbeg, mend, k;
            info.GetLocalRegion(info.GetRank(), mbeg, mend);
            for (k = mbeg; k < mend; ++k) {
                aresid += (test[k] - b[k]) * (test[k] - b[k]);
                bresid += b[k] * b[k];
            }
        }
        double temp[2] = {aresid, bresid}, recv[2] = {aresid, bresid};
Kirill Terekhov's avatar
Kirill Terekhov committed
288
#if defined(USE_MPI)
289
        MPI_Reduce(temp, recv, 2, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD);
Kirill Terekhov's avatar
Kirill Terekhov committed
290
#endif
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
        realresid = sqrt(recv[0] / (recv[1] + 1.0e-100));
        info.RestoreVector(x);
        if (processRank == 0) {
            std::cout << "||Ax-b||=" << sqrt(recv[0]) << " ||b||=" << sqrt(recv[1]) << " ||Ax-b||/||b||=" <<
            sqrt(recv[0] / (recv[1] + 1.0e-100)) << std::endl;
            std::cout << "norms: " << Timer() - tempTimer << std::endl;
            std::cout << processorsCount << " processors for Solver::type=" << type;
            if (success) {
                std::cout << " solved in " << solvingTimer << " secs";
                std::cout << " with " << iters << " iterations to " << resid << " norm";
            } else std::cout << " failed to solve with " << iters << "iterations";
            std::cout << " matrix \"" << matrixFileName << "\"";
            if (vectorBFound) std::cout << " vector \"" << vectorBFileName << "\"";
            std::cout << " true residual ||Ax-b||/||b||=" << realresid;
            std::cout << std::endl;
            std::cout << "reason: " << reason << std::endl;
        }

        if (vectorXFound) {
            Sparse::Vector ex("exact");  // Declare the exact solution vector
            Sparse::Vector err("error"); // Declare the solution error vector
            INMOST_DATA_ENUM_TYPE mbeg, mend, k;
            mat.GetInterval(mbeg, mend);
            err.SetInterval(mbeg, mend);
            ex.Load(vectorXFileName);
            BARRIER
            double dif1 = 0, dif2 = 0, dif8 = 0, norm = 0;
            for (k = mbeg; k < mend; ++k) {
                double dloc = err[k] = fabs(x[k] - ex[k]);
                dif1 += dloc;
                dif2 += dloc * dloc;
                dif8 = (dif8 > dloc) ? dif8 : dloc;
                norm += fabs(ex[k]);
            }
325
#if defined(USE_MPI)
326
327
328
329
330
331
332
333
334
335
            if (processorsCount > 1) {
                double temp1 = dif1;
                MPI_Reduce(&temp1, &dif1, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD);
                temp1 = dif2;
                MPI_Reduce(&temp1, &dif2, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD);
                temp1 = dif8;
                MPI_Reduce(&temp1, &dif8, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD);
                temp1 = norm;
                MPI_Reduce(&temp1, &norm, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD);
            }
336
#endif
337
338
339
340
341
342
343
344
345
346
347
348
349
            dif2 = sqrt(dif2);
            norm += 1.0e-100;
            if (processRank == 0) {
                std::cout << "Difference with exact solution \"" << vectorXFileName << "\": " << std::scientific <<
                std::setprecision(6) << std::endl;
                std::cout << "dif1 = " << dif1 << "  dif2 = " << dif2 << "  dif8 = " << dif8 << "  ||ex||_1 = " <<
                norm << std::endl;
                std::cout << "rel1 = " << dif1 / norm << "  rel2 = " << dif2 / norm << "  rel8 = " << dif8 / norm <<
                std::endl;
            }
        }
    }
    BARRIER
Igor Konshin's avatar
Igor Konshin committed
350
	Solver::Finalize(); // Finalize solver and close MPI activity
Kirill Terekhov's avatar
Kirill Terekhov committed
351
352
	return 0;
}