/[escript]/trunk/paso/profiling/Paso_tests.c
ViewVC logotype

Annotation of /trunk/paso/profiling/Paso_tests.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2160 - (hide annotations)
Tue Dec 16 03:29:54 2008 UTC (14 years, 3 months ago) by artak
File MIME type: text/plain
File size: 4137 byte(s)
Testing envirounment for PASO
1 artak 2160
2     /*******************************************************
3     *
4     * Copyright (c) 2003-2008 by University of Queensland
5     * Earth Systems Science Computational Center (ESSCC)
6     * http://www.uq.edu.au/esscc
7     *
8     * Primary Business: Queensland, Australia
9     * Licensed under the Open Software License version 3.0
10     * http://www.opensource.org/licenses/osl-3.0.php
11     *
12     *******************************************************/
13    
14    
15     /**************************************************************/
16    
17     /* Paso: interface to the direct solvers */
18    
19     /**************************************************************/
20    
21     /* Copyrights by ACcESS Australia 2003 */
22     /* Author: artak@access.edu.au */
23    
24     /**************************************************************/
25    
26     #include "paso/Paso.h"
27     #include "paso/Solver.h"
28     #include "paso/SystemMatrix.h"
29     #include "paso/Options.h"
30     #include "Paso_tests.h"
31    
32    
33     /**************************************************************/
34    
35     void Paso_test_run(Paso_SystemMatrix* A,double* b,dim_t level) {
36    
37     Paso_Options options;
38     Paso_Options_setDefaults(&options);
39    
40     if(level==1) /* Solvers only*/
41     {
42     options.method=PASO_PCG;
43     options.verbose=TRUE;
44     options.preconditioner=PASO_JACOBI;
45     fprintf(stdout,"Test solver: PCG with JACOBI\n");
46     Paso_test_matrix(A,b,&options);
47    
48     fprintf(stdout,"Test solver: GMRES with JACOBI\n");
49     A->solver=NULL;
50     options.method=PASO_GMRES;
51     Paso_test_matrix(A,b,&options);
52    
53     fprintf(stdout,"Test solver: PRES20 with JACOBI\n");
54     A->solver=NULL;
55     options.method=PASO_PRES20;
56     Paso_test_matrix(A,b,&options);
57    
58     fprintf(stdout,"Test solver: MINRES with JACOBI\n");
59     A->solver=NULL;
60     options.method=PASO_MINRES;
61     Paso_test_matrix(A,b,&options);
62    
63     fprintf(stdout,"Test solver: TFQMR with JACOBI\n");
64     A->solver=NULL;
65     options.method=PASO_TFQMR;
66     Paso_test_matrix(A,b,&options);
67     }
68     else if (level==2) /* Preconditiones only with default solver*/
69     {
70     Paso_Options_setDefaults(&options);
71     options.method=PASO_DEFAULT;
72     options.verbose=TRUE;
73     options.preconditioner=PASO_JACOBI;
74     fprintf(stdout,"Test preconditioner: PASO_DEFAULT with JACOBI\n");
75     Paso_test_matrix(A,b,&options);
76    
77     Paso_Options_setDefaults(&options);
78     A->solver=NULL;
79     options.verbose=TRUE;
80     fprintf(stdout,"Test preconditioner: PASO_DEFAULT with ILU\n");
81     options.method=PASO_DEFAULT;
82     options.preconditioner=PASO_ILU0;
83     Paso_test_matrix(A,b,&options);
84    
85     Paso_Options_setDefaults(&options);
86     A->solver=NULL;
87     options.verbose=TRUE;
88     fprintf(stdout,"Test preconditioner: PASO_DEFAULT with RILU\n");
89     options.method=PASO_DEFAULT;
90     options.preconditioner=PASO_RILU;
91     Paso_test_matrix(A,b,&options);
92    
93     Paso_Options_setDefaults(&options);
94     A->solver=NULL;
95     options.verbose=TRUE;
96     fprintf(stdout,"Test preconditioner: PASO_DEFAULT with GS\n");
97     options.method=PASO_DEFAULT;
98     options.preconditioner=PASO_GS;
99     Paso_test_matrix(A,b,&options);
100    
101     Paso_Options_setDefaults(&options);
102     A->solver=NULL;
103     options.verbose=TRUE;
104     fprintf(stdout,"Test preconditioner: PASO_DEFAULT with AMG\n");
105     options.method=PASO_DEFAULT;
106     options.preconditioner=PASO_AMG;
107     Paso_test_matrix(A,b,&options);
108    
109     }
110     }
111    
112     void Paso_test_matrix(Paso_SystemMatrix* A, double* b, Paso_Options* options ) {
113    
114     dim_t n=Paso_SystemMatrix_getTotalNumRows(A);
115     double *out=NULL;
116     out=MEMALLOC(n,double);
117    
118     if (Paso_checkPtr(out)) {
119     fprintf(stderr,"Cannot allocate memory\n");
120     return;
121     }
122     Paso_solve(A,out,b,options);
123    
124     MEMFREE(out);
125    
126     }
127    
128     void Paso_test_data(char *fileName_p, double* b, Paso_Options* options ) {
129    
130     Paso_SystemMatrix* A=NULL;
131     dim_t n=Paso_SystemMatrix_getTotalNumRows(A);
132     double *out=MEMALLOC(n,double);
133     A=Paso_SystemMatrix_loadMM_toCSR(fileName_p);
134    
135     if (Paso_checkPtr(out)) {
136     return;
137     }
138    
139     Paso_solve(A,out,b,options);
140    
141     Paso_SystemMatrix_free(A);
142     MEMFREE(out);
143     }

  ViewVC Help
Powered by ViewVC 1.1.26