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

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

Parent Directory Parent Directory | Revision Log Revision Log


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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26