/[escript]/trunk/paso/src/Solver_ILU.c
ViewVC logotype

Contents of /trunk/paso/src/Solver_ILU.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1841 - (show annotations)
Fri Oct 3 03:57:52 2008 UTC (11 years ago) by gross
File MIME type: text/plain
File size: 25661 byte(s)
modification on LinearPDE class and a first version of Transport class
1
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: ILU preconditioner with reordering */
18
19 /**************************************************************/
20
21 /* Copyrights by ACcESS Australia 2003,2004,2005 */
22 /* Author: gross@access.edu.au */
23
24 /**************************************************************/
25
26 #include "Paso.h"
27 #include "Solver.h"
28 #include "PasoUtil.h"
29
30 /**************************************************************/
31
32 /* free all memory used by ILU */
33
34 void Paso_Solver_ILU_free(Paso_Solver_ILU * in) {
35 if (in!=NULL) {
36 MEMFREE(in->colorOf);
37 MEMFREE(in->factors);
38 MEMFREE(in->main_iptr);
39 Paso_Pattern_free(in->pattern);
40 MEMFREE(in);
41 }
42 }
43
44 /**************************************************************/
45
46 /* constructs the incomplete block factorization of
47
48 */
49 Paso_Solver_ILU* Paso_Solver_getILU(Paso_SparseMatrix * A,bool_t verbose) {
50 dim_t n=A->numRows;
51 dim_t n_block=A->row_block_size;
52 index_t num_colors=0, *mis_marker=NULL;
53 register double A11,A12,A13,A21,A22,A23,A31,A32,A33,D;
54 register double S11,S12,S13,S21,S22,S23,S31,S32,S33;
55 register index_t i,iptr_main,iptr,iptr_ik,k,iptr_kj,j,iptr_ij,color,color2;
56 double time0,time_color,time_fac;
57 /* allocations: */
58 Paso_Solver_ILU* out=MEMALLOC(1,Paso_Solver_ILU);
59 if (Paso_checkPtr(out)) return NULL;
60 out->colorOf=MEMALLOC(n,index_t);
61 out->factors=MEMALLOC(A->len,double);
62 out->main_iptr=MEMALLOC(n,index_t);
63 out->pattern=Paso_Pattern_getReference(A->pattern);
64 out->n_block=n_block;
65 out->n=n;
66 if ( !(Paso_checkPtr(out->colorOf) || Paso_checkPtr(out->main_iptr) || Paso_checkPtr(out->factors)) ) {
67 time0=Paso_timer();
68 Paso_Pattern_color(A->pattern,&out->num_colors,out->colorOf);
69 time_color=Paso_timer()-time0;
70
71 if (Paso_noError()) {
72 time0=Paso_timer();
73 /* find main diagonal and copy matrix values */
74 #pragma omp parallel for schedule(static) private(i,iptr,iptr_main,k)
75 for (i = 0; i < n; ++i) {
76 iptr_main=A->pattern->ptr[0]-1;
77 for (iptr=A->pattern->ptr[i];iptr<A->pattern->ptr[i+1]; iptr++) {
78 if (A->pattern->index[iptr]==i) iptr_main=iptr;
79 for (k=0;k<n_block*n_block;++k) out->factors[n_block*n_block*iptr+k]=A->val[n_block*n_block*iptr+k];
80 }
81 out->main_iptr[i]=iptr_main;
82 if (iptr_main==A->pattern->ptr[0]-1) {
83 Paso_setError(VALUE_ERROR, "Paso_Solver_getILU: no main diagonal");
84 }
85 }
86 /* start factorization */
87
88 #pragma omp barrier
89 for (color=0;color<out->num_colors && Paso_noError();++color) {
90 if (n_block==1) {
91 #pragma omp parallel for schedule(static) private(i,color2,iptr_ik,k,iptr_kj,S11,j,iptr_ij,A11,iptr_main,D)
92 for (i = 0; i < n; ++i) {
93 if (out->colorOf[i]==color) {
94 for (color2=0;color2<color;++color2) {
95 for (iptr_ik=A->pattern->ptr[i];iptr_ik<A->pattern->ptr[i+1]; ++iptr_ik) {
96 k=A->pattern->index[iptr_ik];
97 if (out->colorOf[k]==color2) {
98 A11=out->factors[iptr_ik];
99 /* a_ij=a_ij-a_ik*a_kj */
100 for (iptr_kj=A->pattern->ptr[k];iptr_kj<A->pattern->ptr[k+1]; iptr_kj++) {
101 j=A->pattern->index[iptr_kj];
102 if (out->colorOf[j]>color2) {
103 S11=out->factors[iptr_kj];
104 for (iptr_ij=A->pattern->ptr[i];iptr_ij<A->pattern->ptr[i+1]; iptr_ij++) {
105 if (j==A->pattern->index[iptr_ij]) {
106 out->factors[iptr_ij]-=A11*S11;
107 break;
108 }
109 }
110 }
111 }
112 }
113 }
114 }
115 iptr_main=out->main_iptr[i];
116 D=out->factors[iptr_main];
117 if (ABS(D)>0.) {
118 D=1./D;
119 out->factors[iptr_main]=D;
120 /* a_ik=a_ii^{-1}*a_ik */
121 for (iptr_ik=A->pattern->ptr[i];iptr_ik<A->pattern->ptr[i+1]; ++iptr_ik) {
122 k=A->pattern->index[iptr_ik];
123 if (out->colorOf[k]>color) {
124 A11=out->factors[iptr_ik];
125 out->factors[iptr_ik]=A11*D;
126 }
127 }
128 } else {
129 Paso_setError(ZERO_DIVISION_ERROR, "Paso_Solver_getILU: non-regular main diagonal block.");
130 }
131 }
132 }
133 } else if (n_block==2) {
134 #pragma omp parallel for schedule(static) private(i,color2,iptr_ik,k,iptr_kj,S11,S21,S12,S22,j,iptr_ij,A11,A21,A12,A22,iptr_main,D)
135 for (i = 0; i < n; ++i) {
136 if (out->colorOf[i]==color) {
137 for (color2=0;color2<color;++color2) {
138 for (iptr_ik=A->pattern->ptr[i];iptr_ik<A->pattern->ptr[i+1]; ++iptr_ik) {
139 k=A->pattern->index[iptr_ik];
140 if (out->colorOf[k]==color2) {
141 A11=out->factors[iptr_ik*4 ];
142 A21=out->factors[iptr_ik*4+1];
143 A12=out->factors[iptr_ik*4+2];
144 A22=out->factors[iptr_ik*4+3];
145 /* a_ij=a_ij-a_ik*a_kj */
146 for (iptr_kj=A->pattern->ptr[k];iptr_kj<A->pattern->ptr[k+1]; iptr_kj++) {
147 j=A->pattern->index[iptr_kj];
148 if (out->colorOf[j]>color2) {
149 S11=out->factors[iptr_kj*4];
150 S21=out->factors[iptr_kj*4+1];
151 S12=out->factors[iptr_kj*4+2];
152 S22=out->factors[iptr_kj*4+3];
153 for (iptr_ij=A->pattern->ptr[i];iptr_ij<A->pattern->ptr[i+1]; iptr_ij++) {
154 if (j==A->pattern->index[iptr_ij]) {
155 out->factors[4*iptr_ij ]-=A11*S11+A12*S21;
156 out->factors[4*iptr_ij+1]-=A21*S11+A22*S21;
157 out->factors[4*iptr_ij+2]-=A11*S12+A12*S22;
158 out->factors[4*iptr_ij+3]-=A21*S12+A22*S22;
159 break;
160 }
161 }
162 }
163 }
164 }
165 }
166 }
167 iptr_main=out->main_iptr[i];
168 A11=out->factors[iptr_main*4];
169 A21=out->factors[iptr_main*4+1];
170 A12=out->factors[iptr_main*4+2];
171 A22=out->factors[iptr_main*4+3];
172 D = A11*A22-A12*A21;
173 if (ABS(D)>0.) {
174 D=1./D;
175 S11= A22*D;
176 S21=-A21*D;
177 S12=-A12*D;
178 S22= A11*D;
179 out->factors[iptr_main*4] = S11;
180 out->factors[iptr_main*4+1]= S21;
181 out->factors[iptr_main*4+2]= S12;
182 out->factors[iptr_main*4+3]= S22;
183 /* a_ik=a_ii^{-1}*a_ik */
184 for (iptr_ik=A->pattern->ptr[i];iptr_ik<A->pattern->ptr[i+1]; ++iptr_ik) {
185 k=A->pattern->index[iptr_ik];
186 if (out->colorOf[k]>color) {
187 A11=out->factors[iptr_ik*4 ];
188 A21=out->factors[iptr_ik*4+1];
189 A12=out->factors[iptr_ik*4+2];
190 A22=out->factors[iptr_ik*4+3];
191 out->factors[4*iptr_ik ]=S11*A11+S12*A21;
192 out->factors[4*iptr_ik+1]=S21*A11+S22*A21;
193 out->factors[4*iptr_ik+2]=S11*A12+S12*A22;
194 out->factors[4*iptr_ik+3]=S21*A12+S22*A22;
195 }
196 }
197 } else {
198 Paso_setError(ZERO_DIVISION_ERROR, "Paso_Solver_getILU: non-regular main diagonal block.");
199 }
200 }
201 }
202 } else if (n_block==3) {
203 #pragma omp parallel for schedule(static) private(i,color2,iptr_ik,k,iptr_kj,S11,S21,S31,S12,S22,S32,S13,S23,S33,j,iptr_ij,A11,A21,A31,A12,A22,A32,A13,A23,A33,iptr_main,D)
204 for (i = 0; i < n; ++i) {
205 if (out->colorOf[i]==color) {
206 for (color2=0;color2<color;++color2) {
207 for (iptr_ik=A->pattern->ptr[i];iptr_ik<A->pattern->ptr[i+1]; ++iptr_ik) {
208 k=A->pattern->index[iptr_ik];
209 if (out->colorOf[k]==color2) {
210 A11=out->factors[iptr_ik*9 ];
211 A21=out->factors[iptr_ik*9+1];
212 A31=out->factors[iptr_ik*9+2];
213 A12=out->factors[iptr_ik*9+3];
214 A22=out->factors[iptr_ik*9+4];
215 A32=out->factors[iptr_ik*9+5];
216 A13=out->factors[iptr_ik*9+6];
217 A23=out->factors[iptr_ik*9+7];
218 A33=out->factors[iptr_ik*9+8];
219 /* a_ij=a_ij-a_ik*a_kj */
220 for (iptr_kj=A->pattern->ptr[k];iptr_kj<A->pattern->ptr[k+1]; iptr_kj++) {
221 j=A->pattern->index[iptr_kj];
222 if (out->colorOf[j]>color2) {
223 S11=out->factors[iptr_kj*9 ];
224 S21=out->factors[iptr_kj*9+1];
225 S31=out->factors[iptr_kj*9+2];
226 S12=out->factors[iptr_kj*9+3];
227 S22=out->factors[iptr_kj*9+4];
228 S32=out->factors[iptr_kj*9+5];
229 S13=out->factors[iptr_kj*9+6];
230 S23=out->factors[iptr_kj*9+7];
231 S33=out->factors[iptr_kj*9+8];
232 for (iptr_ij=A->pattern->ptr[i];iptr_ij<A->pattern->ptr[i+1]; iptr_ij++) {
233 if (j==A->pattern->index[iptr_ij]) {
234 out->factors[iptr_ij*9 ]-=A11*S11+A12*S21+A13*S31;
235 out->factors[iptr_ij*9+1]-=A21*S11+A22*S21+A23*S31;
236 out->factors[iptr_ij*9+2]-=A31*S11+A32*S21+A33*S31;
237 out->factors[iptr_ij*9+3]-=A11*S12+A12*S22+A13*S32;
238 out->factors[iptr_ij*9+4]-=A21*S12+A22*S22+A23*S32;
239 out->factors[iptr_ij*9+5]-=A31*S12+A32*S22+A33*S32;
240 out->factors[iptr_ij*9+6]-=A11*S13+A12*S23+A13*S33;
241 out->factors[iptr_ij*9+7]-=A21*S13+A22*S23+A23*S33;
242 out->factors[iptr_ij*9+8]-=A31*S13+A32*S23+A33*S33;
243 break;
244 }
245 }
246 }
247 }
248 }
249 }
250 }
251 iptr_main=out->main_iptr[i];
252 A11=out->factors[iptr_main*9 ];
253 A21=out->factors[iptr_main*9+1];
254 A31=out->factors[iptr_main*9+2];
255 A12=out->factors[iptr_main*9+3];
256 A22=out->factors[iptr_main*9+4];
257 A32=out->factors[iptr_main*9+5];
258 A13=out->factors[iptr_main*9+6];
259 A23=out->factors[iptr_main*9+7];
260 A33=out->factors[iptr_main*9+8];
261 D = A11*(A22*A33-A23*A32)+ A12*(A31*A23-A21*A33)+A13*(A21*A32-A31*A22);
262 if (ABS(D)>0.) {
263 D=1./D;
264 S11=(A22*A33-A23*A32)*D;
265 S21=(A31*A23-A21*A33)*D;
266 S31=(A21*A32-A31*A22)*D;
267 S12=(A13*A32-A12*A33)*D;
268 S22=(A11*A33-A31*A13)*D;
269 S32=(A12*A31-A11*A32)*D;
270 S13=(A12*A23-A13*A22)*D;
271 S23=(A13*A21-A11*A23)*D;
272 S33=(A11*A22-A12*A21)*D;
273
274 out->factors[iptr_main*9 ]=S11;
275 out->factors[iptr_main*9+1]=S21;
276 out->factors[iptr_main*9+2]=S31;
277 out->factors[iptr_main*9+3]=S12;
278 out->factors[iptr_main*9+4]=S22;
279 out->factors[iptr_main*9+5]=S32;
280 out->factors[iptr_main*9+6]=S13;
281 out->factors[iptr_main*9+7]=S23;
282 out->factors[iptr_main*9+8]=S33;
283
284 /* a_ik=a_ii^{-1}*a_ik */
285 for (iptr_ik=A->pattern->ptr[i];iptr_ik<A->pattern->ptr[i+1]; ++iptr_ik) {
286 k=A->pattern->index[iptr_ik];
287 if (out->colorOf[k]>color) {
288 A11=out->factors[iptr_ik*9 ];
289 A21=out->factors[iptr_ik*9+1];
290 A31=out->factors[iptr_ik*9+2];
291 A12=out->factors[iptr_ik*9+3];
292 A22=out->factors[iptr_ik*9+4];
293 A32=out->factors[iptr_ik*9+5];
294 A13=out->factors[iptr_ik*9+6];
295 A23=out->factors[iptr_ik*9+7];
296 A33=out->factors[iptr_ik*9+8];
297 out->factors[iptr_ik*9 ]=S11*A11+S12*A21+S13*A31;
298 out->factors[iptr_ik*9+1]=S21*A11+S22*A21+S23*A31;
299 out->factors[iptr_ik*9+2]=S31*A11+S32*A21+S33*A31;
300 out->factors[iptr_ik*9+3]=S11*A12+S12*A22+S13*A32;
301 out->factors[iptr_ik*9+4]=S21*A12+S22*A22+S23*A32;
302 out->factors[iptr_ik*9+5]=S31*A12+S32*A22+S33*A32;
303 out->factors[iptr_ik*9+6]=S11*A13+S12*A23+S13*A33;
304 out->factors[iptr_ik*9+7]=S21*A13+S22*A23+S23*A33;
305 out->factors[iptr_ik*9+8]=S31*A13+S32*A23+S33*A33;
306 }
307 }
308 } else {
309 Paso_setError(ZERO_DIVISION_ERROR, "Paso_Solver_getILU: non-regular main diagonal block.");
310 }
311 }
312 }
313 } else {
314 Paso_setError(VALUE_ERROR, "Paso_Solver_getILU: block size greater than 3 is not supported.");
315 }
316 #pragma omp barrier
317 }
318 time_fac=Paso_timer()-time0;
319 }
320 }
321 if (Paso_noError()) {
322 if (verbose) {
323 printf("ILU: %d color used \n",out->num_colors);
324 printf("timing: ILU: coloring/elemination : %e/%e\n",time_color,time_fac);
325 }
326 return out;
327 } else {
328 Paso_Solver_ILU_free(out);
329 return NULL;
330 }
331 }
332
333 /**************************************************************/
334
335 /* apply ILU precondition b-> x
336
337 in fact it solves LUx=b in the form x= U^{-1} L^{-1}b
338
339 should be called within a parallel region
340 barrier synconization should be performed to make sure that the input vector available
341
342 */
343
344 void Paso_Solver_solveILU(Paso_Solver_ILU * ilu, double * x, double * b) {
345 register dim_t i,k;
346 register index_t color,iptr_ik,iptr_main;
347 register double S1,S2,S3,R1,R2,R3;
348 dim_t n_block=ilu->n_block;
349 dim_t n=ilu->n;
350
351
352 /* copy x into b*/
353 #pragma omp parallel for private(i) schedule(static)
354 for (i=0;i<n*n_block;++i) x[i]=b[i];
355 /* forward substitution */
356 for (color=0;color<ilu->num_colors;++color) {
357 if (n_block==1) {
358 #pragma omp parallel for schedule(static) private(i,iptr_ik,k,S1,R1,iptr_main)
359 for (i = 0; i < n; ++i) {
360 if (ilu->colorOf[i]==color) {
361 /* x_i=x_i-a_ik*x_k */
362 S1=x[i];
363 for (iptr_ik=ilu->pattern->ptr[i];iptr_ik<ilu->pattern->ptr[i+1]; ++iptr_ik) {
364 k=ilu->pattern->index[iptr_ik];
365 if (ilu->colorOf[k]<color) {
366 R1=x[k];
367 S1-=ilu->factors[iptr_ik]*R1;
368 }
369 }
370 iptr_main=ilu->main_iptr[i];
371 x[i]=ilu->factors[iptr_main]*S1;
372 }
373 }
374 } else if (n_block==2) {
375 #pragma omp parallel for schedule(static) private(i,iptr_ik,k,iptr_main,S1,S2,R1,R2)
376 for (i = 0; i < n; ++i) {
377 if (ilu->colorOf[i]==color) {
378 /* x_i=x_i-a_ik*x_k */
379 S1=x[2*i];
380 S2=x[2*i+1];
381 for (iptr_ik=ilu->pattern->ptr[i];iptr_ik<ilu->pattern->ptr[i+1]; ++iptr_ik) {
382 k=ilu->pattern->index[iptr_ik];
383 if (ilu->colorOf[k]<color) {
384 R1=x[2*k];
385 R2=x[2*k+1];
386 S1-=ilu->factors[4*iptr_ik ]*R1+ilu->factors[4*iptr_ik+2]*R2;
387 S2-=ilu->factors[4*iptr_ik+1]*R1+ilu->factors[4*iptr_ik+3]*R2;
388 }
389 }
390 iptr_main=ilu->main_iptr[i];
391 x[2*i ]=ilu->factors[4*iptr_main ]*S1+ilu->factors[4*iptr_main+2]*S2;
392 x[2*i+1]=ilu->factors[4*iptr_main+1]*S1+ilu->factors[4*iptr_main+3]*S2;
393 }
394
395 }
396 } else if (n_block==3) {
397 #pragma omp parallel for schedule(static) private(i,iptr_ik,iptr_main,k,S1,S2,S3,R1,R2,R3)
398 for (i = 0; i < n; ++i) {
399 if (ilu->colorOf[i]==color) {
400 /* x_i=x_i-a_ik*x_k */
401 S1=x[3*i];
402 S2=x[3*i+1];
403 S3=x[3*i+2];
404 for (iptr_ik=ilu->pattern->ptr[i];iptr_ik<ilu->pattern->ptr[i+1]; ++iptr_ik) {
405 k=ilu->pattern->index[iptr_ik];
406 if (ilu->colorOf[k]<color) {
407 R1=x[3*k];
408 R2=x[3*k+1];
409 R3=x[3*k+2];
410 S1-=ilu->factors[9*iptr_ik ]*R1+ilu->factors[9*iptr_ik+3]*R2+ilu->factors[9*iptr_ik+6]*R3;
411 S2-=ilu->factors[9*iptr_ik+1]*R1+ilu->factors[9*iptr_ik+4]*R2+ilu->factors[9*iptr_ik+7]*R3;
412 S3-=ilu->factors[9*iptr_ik+2]*R1+ilu->factors[9*iptr_ik+5]*R2+ilu->factors[9*iptr_ik+8]*R3;
413 }
414 }
415 iptr_main=ilu->main_iptr[i];
416 x[3*i ]=ilu->factors[9*iptr_main ]*S1+ilu->factors[9*iptr_main+3]*S2+ilu->factors[9*iptr_main+6]*S3;
417 x[3*i+1]=ilu->factors[9*iptr_main+1]*S1+ilu->factors[9*iptr_main+4]*S2+ilu->factors[9*iptr_main+7]*S3;
418 x[3*i+2]=ilu->factors[9*iptr_main+2]*S1+ilu->factors[9*iptr_main+5]*S2+ilu->factors[9*iptr_main+8]*S3;
419 }
420 }
421 }
422 }
423 /* backward substitution */
424 for (color=(ilu->num_colors)-1;color>-1;--color) {
425 if (n_block==1) {
426 #pragma omp parallel for schedule(static) private(i,iptr_ik,k,S1,R1)
427 for (i = 0; i < n; ++i) {
428 if (ilu->colorOf[i]==color) {
429 /* x_i=x_i-a_ik*x_k */
430 S1=x[i];
431 for (iptr_ik=ilu->pattern->ptr[i];iptr_ik<ilu->pattern->ptr[i+1]; ++iptr_ik) {
432 k=ilu->pattern->index[iptr_ik];
433 if (ilu->colorOf[k]>color) {
434 R1=x[k];
435 S1-=ilu->factors[iptr_ik]*R1;
436 }
437 }
438 x[i]=S1;
439 }
440 }
441 } else if (n_block==2) {
442 #pragma omp parallel for schedule(static) private(i,iptr_ik,k,S1,S2,R1,R2)
443 for (i = 0; i < n; ++i) {
444 if (ilu->colorOf[i]==color) {
445 /* x_i=x_i-a_ik*x_k */
446 S1=x[2*i];
447 S2=x[2*i+1];
448 for (iptr_ik=ilu->pattern->ptr[i];iptr_ik<ilu->pattern->ptr[i+1]; ++iptr_ik) {
449 k=ilu->pattern->index[iptr_ik];
450 if (ilu->colorOf[k]>color) {
451 R1=x[2*k];
452 R2=x[2*k+1];
453 S1-=ilu->factors[4*iptr_ik ]*R1+ilu->factors[4*iptr_ik+2]*R2;
454 S2-=ilu->factors[4*iptr_ik+1]*R1+ilu->factors[4*iptr_ik+3]*R2;
455 }
456 }
457 x[2*i]=S1;
458 x[2*i+1]=S2;
459 }
460 }
461 } else if (n_block==3) {
462 #pragma omp parallel for schedule(static) private(i,iptr_ik,k,S1,S2,S3,R1,R2,R3)
463 for (i = 0; i < n; ++i) {
464 if (ilu->colorOf[i]==color) {
465 /* x_i=x_i-a_ik*x_k */
466 S1=x[3*i ];
467 S2=x[3*i+1];
468 S3=x[3*i+2];
469 for (iptr_ik=ilu->pattern->ptr[i];iptr_ik<ilu->pattern->ptr[i+1]; ++iptr_ik) {
470 k=ilu->pattern->index[iptr_ik];
471 if (ilu->colorOf[k]>color) {
472 R1=x[3*k];
473 R2=x[3*k+1];
474 R3=x[3*k+2];
475 S1-=ilu->factors[9*iptr_ik ]*R1+ilu->factors[9*iptr_ik+3]*R2+ilu->factors[9*iptr_ik+6]*R3;
476 S2-=ilu->factors[9*iptr_ik+1]*R1+ilu->factors[9*iptr_ik+4]*R2+ilu->factors[9*iptr_ik+7]*R3;
477 S3-=ilu->factors[9*iptr_ik+2]*R1+ilu->factors[9*iptr_ik+5]*R2+ilu->factors[9*iptr_ik+8]*R3;
478 }
479 }
480 x[3*i]=S1;
481 x[3*i+1]=S2;
482 x[3*i+2]=S3;
483 }
484 }
485 }
486 #pragma omp barrier
487 }
488 return;
489 }
490 /*
491 * $Log$
492 * Revision 1.2 2005/09/15 03:44:40 jgs
493 * Merge of development branch dev-02 back to main trunk on 2005-09-15
494 *
495 * Revision 1.1.2.1 2005/09/05 06:29:50 gross
496 * These files have been extracted from finley to define a stand alone libray for iterative
497 * linear solvers on the ALTIX. main entry through Paso_solve. this version compiles but
498 * has not been tested yet.
499 *
500 *
501 */

Properties

Name Value
svn:eol-style native
svn:keywords Author Date Id Revision

  ViewVC Help
Powered by ViewVC 1.1.26