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

Diff of /trunk/paso/src/Solver_AMG.c

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

revision 1890 by artak, Fri Oct 17 00:14:22 2008 UTC revision 2551 by gross, Thu Jul 23 09:19:15 2009 UTC
# Line 1  Line 1 
1    
2  /*******************************************************  /*******************************************************
3  *  *
4  * Copyright (c) 2003-2008 by University of Queensland  * Copyright (c) 2003-2009 by University of Queensland
5  * Earth Systems Science Computational Center (ESSCC)  * Earth Systems Science Computational Center (ESSCC)
6  * http://www.uq.edu.au/esscc  * http://www.uq.edu.au/esscc
7  *  *
# Line 14  Line 14 
14    
15  /**************************************************************/  /**************************************************************/
16    
17  /* Paso: AMG preconditioner with reordering                 */  /* Paso: AMG preconditioner                                  */
18    
19  /**************************************************************/  /**************************************************************/
20    
# Line 24  Line 24 
24    
25  #include "Paso.h"  #include "Paso.h"
26  #include "Solver.h"  #include "Solver.h"
27    #include "Options.h"
28  #include "PasoUtil.h"  #include "PasoUtil.h"
29    #include "UMFPACK.h"
30    #include "MKL.h"
31    #include "SystemMatrix.h"
32    #include "Pattern_coupling.h"
33    
34  /**************************************************************/  /**************************************************************/
35    
# Line 33  Line 38 
38  void Paso_Solver_AMG_free(Paso_Solver_AMG * in) {  void Paso_Solver_AMG_free(Paso_Solver_AMG * in) {
39       if (in!=NULL) {       if (in!=NULL) {
40          Paso_Solver_AMG_free(in->AMG_of_Schur);          Paso_Solver_AMG_free(in->AMG_of_Schur);
41            Paso_Solver_Jacobi_free(in->GS);
42          MEMFREE(in->inv_A_FF);          MEMFREE(in->inv_A_FF);
43          MEMFREE(in->A_FF_pivot);          MEMFREE(in->A_FF_pivot);
44          Paso_SparseMatrix_free(in->A_FC);          Paso_SparseMatrix_free(in->A_FC);
# Line 45  void Paso_Solver_AMG_free(Paso_Solver_AM Line 51  void Paso_Solver_AMG_free(Paso_Solver_AM
51          MEMFREE(in->b_F);          MEMFREE(in->b_F);
52          MEMFREE(in->x_C);          MEMFREE(in->x_C);
53          MEMFREE(in->b_C);          MEMFREE(in->b_C);
54            #ifdef MKL
55            #else
56               #ifdef UMFPACK
57               Paso_UMFPACK1_free((Paso_UMFPACK_Handler*)(in->solver));
58               #endif
59            #endif
60            in->solver=NULL;
61            MEMFREE(in->b_C);
62          MEMFREE(in);          MEMFREE(in);
63       }       }
64  }  }
# Line 68  to Line 82  to
82     then AMG is applied to S again until S becomes empty     then AMG is applied to S again until S becomes empty
83    
84  */  */
85  Paso_Solver_AMG* Paso_Solver_getAMG(Paso_SparseMatrix *A_p,bool_t verbose,dim_t level) {  Paso_Solver_AMG* Paso_Solver_getAMG(Paso_SparseMatrix *A_p,dim_t level,Paso_Options* options) {
86    Paso_Solver_AMG* out=NULL;    Paso_Solver_AMG* out=NULL;
87      bool_t verbose=options->verbose;
88    dim_t n=A_p->numRows;    dim_t n=A_p->numRows;
89    dim_t n_block=A_p->row_block_size;    dim_t n_block=A_p->row_block_size;
90    index_t* mis_marker=NULL;      index_t* mis_marker=NULL;  
91    index_t* counter=NULL;    index_t* counter=NULL;
   double *rs=NULL;    
92    index_t iPtr,*index, *where_p, iPtr_s;    index_t iPtr,*index, *where_p, iPtr_s;
93    dim_t i,k,j,j0;    dim_t i,k,j;
94    Paso_SparseMatrix * schur=NULL;    Paso_SparseMatrix * schur=NULL;
95    Paso_SparseMatrix * schur_withFillIn=NULL;    Paso_SparseMatrix * schur_withFillIn=NULL;
96    schur_withFillIn=MEMALLOC(1,Paso_SparseMatrix);    double S=0;
     
97        
98    double A11,A12,A13,A21,A22,A23,A31,A32,A33,D,time0,time1,time2;    /*Make sure we have block sizes 1*/
99        if (A_p->col_block_size>1) {
100         Paso_setError(TYPE_ERROR,"Paso_Solver_getAMG: AMG requires column block size 1.");
101         return NULL;
102      }
103      if (A_p->row_block_size>1) {
104         Paso_setError(TYPE_ERROR,"Paso_Solver_getAMG: AMG requires row block size 1.");
105         return NULL;
106      }
107      out=MEMALLOC(1,Paso_Solver_AMG);
108    /* identify independend set of rows/columns */    /* identify independend set of rows/columns */
109    mis_marker=TMPMEMALLOC(n,index_t);    mis_marker=TMPMEMALLOC(n,index_t);
110    counter=TMPMEMALLOC(n,index_t);    counter=TMPMEMALLOC(n,index_t);
111    rs=TMPMEMALLOC(n,double);    if ( !( Paso_checkPtr(mis_marker) || Paso_checkPtr(counter) || Paso_checkPtr(out)) ) {
112    out=MEMALLOC(1,Paso_Solver_AMG);       out->AMG_of_Schur=NULL;
113    out->AMG_of_Schur=NULL;       out->inv_A_FF=NULL;
114    out->inv_A_FF=NULL;       out->A_FF_pivot=NULL;
115    out->A_FF_pivot=NULL;       out->A_FC=NULL;
116    out->A_FC=NULL;       out->A_CF=NULL;
117    out->A_CF=NULL;       out->rows_in_F=NULL;
118    out->rows_in_F=NULL;       out->rows_in_C=NULL;
119    out->rows_in_C=NULL;       out->mask_F=NULL;
120    out->mask_F=NULL;       out->mask_C=NULL;
121    out->mask_C=NULL;       out->x_F=NULL;
122    out->x_F=NULL;       out->b_F=NULL;
123    out->b_F=NULL;       out->x_C=NULL;
124    out->x_C=NULL;       out->b_C=NULL;
125    out->b_C=NULL;       out->GS=NULL;
126    out->A=Paso_SparseMatrix_getReference(A_p);       out->A=Paso_SparseMatrix_getReference(A_p);
127    out->level=level;       out->GS=NULL;
128         out->solver=NULL;
129         /*out->GS=Paso_Solver_getGS(A_p,verbose);*/
130         out->level=level;
131        
132   /* fprintf(stderr,"START OF MATRIX \n\n");       if (level==0 || n<=options->min_coarse_matrix_size) {
133    for (i = 0; i < A_p->numRows; ++i) {           out->coarsest_level=TRUE;
134       for (iPtr=A_p->pattern->ptr[i];iPtr<A_p->pattern->ptr[i + 1]; ++iPtr) {           #ifdef MKL
135         j=A_p->pattern->index[iPtr];           #else
136         fprintf(stderr,"A[%d,%d]=%.2f ",i,j,A_p->val[iPtr]);              #ifdef UMFPACK
137       }              #else
138       fprintf(stderr,"\n");                  out->GS=Paso_Solver_getJacobi(A_p);
139     }              #endif
140     fprintf(stderr,"END OF MATRIX \n\n");           #endif
141   */       } else {
142    if ( !(Paso_checkPtr(mis_marker) || Paso_checkPtr(out) || Paso_checkPtr(counter) ) ) {           out->coarsest_level=FALSE;
143       /* identify independend set of rows/columns */           out->GS=Paso_Solver_getJacobi(A_p);
144       time0=Paso_timer();  
145       #pragma omp parallel for private(i) schedule(static)           /* identify independend set of rows/columns */
146       for (i=0;i<n;++i) mis_marker[i]=-1;           #pragma omp parallel for private(i) schedule(static)
147       Paso_Pattern_RS(A_p,mis_marker,0.25);           for (i=0;i<n;++i) mis_marker[i]=-1;
148       /*  
149       for (i=0;i<n;++i) fprintf(stderr," i=%d mis[i]=%d \n",i,mis_marker[i]);           if (options->coarsening_method == PASO_YAIR_SHAPIRA_COARSENING) {
150       */                Paso_Pattern_coup(A_p,mis_marker,options->coarsening_threshold);
151       time2=Paso_timer()-time0;           }
152       if (Paso_noError()) {           else if (options->coarsening_method == PASO_RUGE_STUEBEN_COARSENING) {
153          #pragma omp parallel for private(i) schedule(static)                Paso_Pattern_RS(A_p,mis_marker,options->coarsening_threshold);
154          for (i = 0; i < n; ++i) counter[i]=mis_marker[i];           }
155          out->n=n;           else if (options->coarsening_method == PASO_AGGREGATION_COARSENING) {
156          out->n_block=n_block;               Paso_Pattern_Aggregiation(A_p,mis_marker,options->coarsening_threshold);
157          out->n_F=Paso_Util_cumsum(n,counter);          }
158          out->mask_F=MEMALLOC(n,index_t);          else {
159          out->rows_in_F=MEMALLOC(out->n_F,index_t);             /*Default coarseneing*/
160          out->inv_A_FF=MEMALLOC(n_block*n_block*out->n_F,double);             /* Paso_Pattern_RS(A_p,mis_marker,options->coarsening_threshold); */
161          out->A_FF_pivot=NULL; /* later use for block size>3 */  
162          if (! (Paso_checkPtr(out->mask_F) || Paso_checkPtr(out->inv_A_FF) || Paso_checkPtr(out->rows_in_F) ) ) {               Paso_Pattern_Aggregiation(A_p,mis_marker,options->coarsening_threshold);
163             #pragma omp parallel          }
164             {      
165            if (Paso_noError()) {
166               #pragma omp parallel for private(i) schedule(static)
167               for (i = 0; i < n; ++i) counter[i]=mis_marker[i];
168    
169               out->n=n;
170               out->n_block=n_block;
171               out->n_F=Paso_Util_cumsum(n,counter);
172               out->mask_F=MEMALLOC(n,index_t);
173               out->rows_in_F=MEMALLOC(out->n_F,index_t);
174               out->inv_A_FF=MEMALLOC(n_block*n_block*out->n_F,double);
175               out->A_FF_pivot=NULL; /* later use for block size>3 */
176               if (! (Paso_checkPtr(out->mask_F) || Paso_checkPtr(out->inv_A_FF) || Paso_checkPtr(out->rows_in_F) ) ) {
177                /* creates an index for F from mask */                /* creates an index for F from mask */
178                #pragma omp for private(i) schedule(static)                #pragma omp parallel for private(i) schedule(static)
179                for (i = 0; i < out->n_F; ++i) out->rows_in_F[i]=-1;                for (i = 0; i < out->n_F; ++i) out->rows_in_F[i]=-1;
180                #pragma omp for private(i) schedule(static)                #pragma omp parallel for private(i) schedule(static)
181                for (i = 0; i < n; ++i) {                for (i = 0; i < n; ++i) {
182                   if  (mis_marker[i]) {                   if  (mis_marker[i]) {
183                          out->rows_in_F[counter[i]]=i;                          out->rows_in_F[counter[i]]=i;
# Line 151  Paso_Solver_AMG* Paso_Solver_getAMG(Paso Line 186  Paso_Solver_AMG* Paso_Solver_getAMG(Paso
186                          out->mask_F[i]=-1;                          out->mask_F[i]=-1;
187                   }                   }
188                }                }
189                /* Compute row-sum for getting rs(A_FF)*/                
190                #pragma omp for private(i,iPtr) schedule(static)                /* Compute row-sum for getting rs(A_FF)^-1*/
191                  #pragma omp parallel for private(i,iPtr,j,S) schedule(static)
192                for (i = 0; i < out->n_F; ++i) {                for (i = 0; i < out->n_F; ++i) {
193                  rs[i]=0;                  S=0;
194    printf("%d : %d -> ",i, out->rows_in_F[i]);
195                  for (iPtr=A_p->pattern->ptr[out->rows_in_F[i]];iPtr<A_p->pattern->ptr[out->rows_in_F[i] + 1]; ++iPtr) {                  for (iPtr=A_p->pattern->ptr[out->rows_in_F[i]];iPtr<A_p->pattern->ptr[out->rows_in_F[i] + 1]; ++iPtr) {
196                   rs[i]+=A_p->val[iPtr];                   j=A_p->pattern->index[iPtr];
197    if (j==out->rows_in_F[i]) printf("%e",A_p->val[iPtr]);
198                     if (mis_marker[j])
199                         S+=A_p->val[iPtr];
200                  }                  }
201    printf("-> %e\n",S);
202                    out->inv_A_FF[i]=1./S;
203                }                }
204                           }
205                #pragma omp for private(i, where_p,iPtr,A11,A12,A13,A21,A22,A23,A31,A32,A33,D,index) schedule(static)          }
               for (i = 0; i < out->n_F; i++) {  
                 /* find main diagonal */  
                 iPtr=A_p->pattern->ptr[out->rows_in_F[i]];  
                 index=&(A_p->pattern->index[iPtr]);  
                 where_p=(index_t*)bsearch(&out->rows_in_F[i],  
                                         index,  
                                         A_p->pattern->ptr[out->rows_in_F[i] + 1]-A_p->pattern->ptr[out->rows_in_F[i]],  
                                         sizeof(index_t),  
                                         Paso_comparIndex);  
                 if (where_p==NULL) {  
                     Paso_setError(VALUE_ERROR, "Paso_Solver_getAMG: main diagonal element missing.");  
                 } else {  
                     iPtr+=(index_t)(where_p-index);  
                     /* get inverse of A_FF block: */  
                       if (ABS(rs[i])>0.) {  
                             out->inv_A_FF[i]=1./rs[i];  
                       } else {  
                         out->inv_A_FF[i]=0;  
                       }  
                       
                      /* } else {  
                             Paso_setError(ZERO_DIVISION_ERROR, "Paso_Solver_getAMG: Break-down in AMG decomposition: non-regular main diagonal block.");  
                       }*/  
                 }  
               }  
            } /* end parallel region */  
206    
207             if( Paso_noError()) {          if ( Paso_noError()) {
208                /* if there are no nodes in the coarse level there is no more work to do */                /* if there are no nodes in the coarse level there is no more work to do */
209                out->n_C=n-out->n_F;                out->n_C=n-out->n_F;
210                 /*if (out->n_C>11) {*/              
211                 if (level>0) {               /*if (out->n_F>500) */
212                     out->rows_in_C=MEMALLOC(out->n_C,index_t);                out->rows_in_C=MEMALLOC(out->n_C,index_t);
213                     out->mask_C=MEMALLOC(n,index_t);                out->mask_C=MEMALLOC(n,index_t);
214                     if (! (Paso_checkPtr(out->mask_C) || Paso_checkPtr(out->rows_in_C) ) ) {                if (! (Paso_checkPtr(out->mask_C) || Paso_checkPtr(out->rows_in_C) ) ) {
215                         /* creates an index for C from mask */                     /* creates an index for C from mask */
216                         #pragma omp parallel for private(i) schedule(static)                     #pragma omp parallel for private(i) schedule(static)
217                         for (i = 0; i < n; ++i) counter[i]=! mis_marker[i];                     for (i = 0; i < n; ++i) counter[i]=! mis_marker[i];
218                         Paso_Util_cumsum(n,counter);                     Paso_Util_cumsum(n,counter);
219                         #pragma omp parallel                     #pragma omp parallel for private(i) schedule(static)
220                         {                     for (i = 0; i < out->n_C; ++i) out->rows_in_C[i]=-1;
221                            #pragma omp for private(i) schedule(static)                     #pragma omp parallel for private(i) schedule(static)
222                            for (i = 0; i < out->n_C; ++i) out->rows_in_C[i]=-1;                     for (i = 0; i < n; ++i) {
223                            #pragma omp for private(i) schedule(static)                              if  (! mis_marker[i]) {
                           for (i = 0; i < n; ++i) {  
                              if  (! mis_marker[i]) {  
224                                  out->rows_in_C[counter[i]]=i;                                  out->rows_in_C[counter[i]]=i;
225                                  out->mask_C[i]=counter[i];                                  out->mask_C[i]=counter[i];
226                               } else {                               } else {
227                                  out->mask_C[i]=-1;                                  out->mask_C[i]=-1;
228                               }                               }
229                            }                     }
230                        } /* end parallel region */                }
231                        /* get A_CF block: */          }
232                        out->A_CF=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_F,out->rows_in_C,out->mask_F);          if ( Paso_noError()) {
233                        if (Paso_noError()) {                  /* get A_CF block: */
234                           /* get A_FC block: */                  out->A_CF=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_F,out->rows_in_C,out->mask_F);
235                           out->A_FC=Paso_SparseMatrix_getSubmatrix(A_p,out->n_F,out->n_C,out->rows_in_F,out->mask_C);                  /* get A_FC block: */
236                           /* get A_CC block: */                  out->A_FC=Paso_SparseMatrix_getSubmatrix(A_p,out->n_F,out->n_C,out->rows_in_F,out->mask_C);
237                           if (Paso_noError()) {                  /* get A_CC block: */
238                              schur=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_C,out->rows_in_C,out->mask_C);                  schur=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_C,out->rows_in_C,out->mask_C);
239                                
240                              /*find the pattern of the schur complement with fill in*/          }
241                              schur_withFillIn=Paso_SparseMatrix_alloc(A_p->type,Paso_Pattern_binop(PATTERN_FORMAT_DEFAULT, schur->pattern, Paso_Pattern_multiply(PATTERN_FORMAT_DEFAULT,out->A_CF->pattern,out->A_FC->pattern)),1,1);          if ( Paso_noError()) {
242                                               /*find the pattern of the schur complement with fill in*/
243                              /* copy values over*/      
244                              #pragma omp for private(i,iPtr,iPtr_s,j,j0) schedule(static)                schur_withFillIn=Paso_SparseMatrix_alloc(A_p->type,Paso_Pattern_binop(PATTERN_FORMAT_DEFAULT, schur->pattern, Paso_Pattern_multiply(PATTERN_FORMAT_DEFAULT,out->A_CF->pattern,out->A_FC->pattern)),1,1, TRUE);
245                              for (i = 0; i < schur_withFillIn->numRows; ++i) {          }
246                                for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {          if ( Paso_noError()) {
247                                  j=schur_withFillIn->pattern->index[iPtr];                /* copy values over*/
248                                  schur_withFillIn->val[iPtr]=0.;                #pragma omp parallel for private(i,iPtr,j,iPtr_s,index,where_p) schedule(static)
249                                  for (iPtr_s=schur->pattern->ptr[i];iPtr_s<schur->pattern->ptr[i + 1]; ++iPtr_s){                for (i = 0; i < schur_withFillIn->numRows; ++i) {
250                                      j0=schur->pattern->index[iPtr_s];                  for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {
251                                      if (j==j0) {                     j=schur_withFillIn->pattern->index[iPtr];
252                                        schur_withFillIn->val[iPtr]=schur->val[iPtr_s];                     iPtr_s=schur->pattern->ptr[i];
253                                        break;                     schur_withFillIn->val[iPtr]=0.;
254                                      }                     index=&(schur->pattern->index[iPtr_s]);
255                                  }                     where_p=(index_t*)bsearch(&j,
256                                }                                          index,
257                              }                                          schur->pattern->ptr[i + 1]-schur->pattern->ptr[i],
258                                                                        sizeof(index_t),
259                            /*  for (i = 0; i < schur_withFillIn->numRows; ++i) {                                          Paso_comparIndex);
260                                for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {                     if (where_p!=NULL) {
261                                  j=schur_withFillIn->pattern->index[iPtr];                            schur_withFillIn->val[iPtr]=schur->val[iPtr_s+(index_t)(where_p-index)];
262                                  fprintf(stderr,"A_CC[%d,%d]=%.2f ",i,j,schur_withFillIn->val[iPtr]);                     }
                               }  
                               fprintf(stderr,"\n");  
                             }*/  
                             time0=Paso_timer()-time0;  
                             if (Paso_noError()) {  
                                 time1=Paso_timer();  
                                 /* update A_CC block to get Schur complement and then apply AMG to it */  
                                 Paso_Solver_updateIncompleteSchurComplement(schur_withFillIn,out->A_CF,out->inv_A_FF,out->A_FF_pivot,out->A_FC);  
                                 time1=Paso_timer()-time1;  
                                 out->AMG_of_Schur=Paso_Solver_getAMG(schur_withFillIn,verbose,level-1);  
                                   
                                 Paso_SparseMatrix_free(schur);  
                                 /* Paso_SparseMatrix_free(schur_withFillIn);*/  
                             }  
                             /* allocate work arrays for AMG application */  
                             if (Paso_noError()) {  
                               out->x_F=MEMALLOC(n_block*out->n_F,double);  
                               out->b_F=MEMALLOC(n_block*out->n_F,double);  
                               out->x_C=MEMALLOC(n_block*out->n_C,double);  
                               out->b_C=MEMALLOC(n_block*out->n_C,double);  
   
                               if (! (Paso_checkPtr(out->x_F) || Paso_checkPtr(out->b_F) || Paso_checkPtr(out->x_C) || Paso_checkPtr(out->b_C) ) ) {  
                                   #pragma omp parallel  
                                   {  
                                     #pragma omp for private(i,k) schedule(static)  
                                     for (i = 0; i < out->n_F; ++i) {  
                                           for (k=0; k<n_block;++k) {  
                                              out->x_F[i*n_block+k]=0.;  
                                              out->b_F[i*n_block+k]=0.;  
                                           }  
                                     }  
                                     #pragma omp for private(i,k) schedule(static)  
                                     for (i = 0; i < out->n_C; ++i) {  
                                         for (k=0; k<n_block;++k) {  
                                           out->x_C[i*n_block+k]=0.;  
                                           out->b_C[i*n_block+k]=0.;  
                                         }  
                                     }  
                                   } /* end parallel region */  
                               }  
                             }  
                          }  
                      }  
263                   }                   }
264                }                }
265             }                Paso_Solver_updateIncompleteSchurComplement(schur_withFillIn,out->A_CF,out->inv_A_FF,out->A_FF_pivot,out->A_FC);
266                  out->AMG_of_Schur=Paso_Solver_getAMG(schur_withFillIn,level-1,options);
267            }
268            /* allocate work arrays for AMG application */
269            if (Paso_noError()) {
270                       out->x_F=MEMALLOC(n_block*out->n_F,double);
271                       out->b_F=MEMALLOC(n_block*out->n_F,double);
272                       out->x_C=MEMALLOC(n_block*out->n_C,double);
273                       out->b_C=MEMALLOC(n_block*out->n_C,double);
274    
275                       if (! (Paso_checkPtr(out->x_F) || Paso_checkPtr(out->b_F) || Paso_checkPtr(out->x_C) || Paso_checkPtr(out->b_C) ) ) {
276                           #pragma omp parallel for private(i,k) schedule(static)
277                           for (i = 0; i < out->n_F; ++i) {
278                                for (k=0; k<n_block;++k) {
279                                       out->x_F[i*n_block+k]=0.;
280                                       out->b_F[i*n_block+k]=0.;
281                                }
282                            }
283                            #pragma omp parallel for private(i,k) schedule(static)
284                            for (i = 0; i < out->n_C; ++i) {
285                               for (k=0; k<n_block;++k) {
286                                   out->x_C[i*n_block+k]=0.;
287                                   out->b_C[i*n_block+k]=0.;
288                               }
289                            }
290                       }
291          }          }
292       }       }
293    }    }
294    TMPMEMFREE(mis_marker);    TMPMEMFREE(mis_marker);
295    TMPMEMFREE(counter);    TMPMEMFREE(counter);
296    TMPMEMFREE(rs);    Paso_SparseMatrix_free(schur);
297      Paso_SparseMatrix_free(schur_withFillIn);
298    if (Paso_noError()) {    if (Paso_noError()) {
299        if (verbose) {        if (verbose && level>0 && !out->coarsest_level) {
300           printf("AMG: %d unknowns eliminated. %d left.\n",out->n_F,n-out->n_F);           printf("AMG: level: %d: %d unknowns eliminated. %d left.\n",level, out->n_F,out->n_C);
          if (out->n_C>0) {  
             printf("timing: AMG: MIS/reordering/elemination : %e/%e/%e\n",time2,time0,time1);  
          } else {  
             printf("timing: AMG: MIS: %e\n",time2);  
          }  
301       }       }
302       return out;       return out;
303    } else  {    } else  {
# Line 320  Paso_Solver_AMG* Paso_Solver_getAMG(Paso Line 312  Paso_Solver_AMG* Paso_Solver_getAMG(Paso
312    
313       in fact it solves       in fact it solves
314    
315    [      I         0  ]  [ A_FF 0 ] [ I    invA_FF*A_FF ]  [ x_F ]  = [b_F]    [      I         0  ]  [ A_FF 0 ] [ I    invA_FF*A_FC ]  [ x_F ]  = [b_F]
316    [ A_CF*invA_FF   I  ]  [   0  S ] [ 0          I      ]  [ x_C ]  = [b_C]    [ A_CF*invA_FF   I  ]  [   0  S ] [ 0          I      ]  [ x_C ]  = [b_C]
317    
318   in the form   in the form
# Line 339  Paso_Solver_AMG* Paso_Solver_getAMG(Paso Line 331  Paso_Solver_AMG* Paso_Solver_getAMG(Paso
331  */  */
332    
333  void Paso_Solver_solveAMG(Paso_Solver_AMG * amg, double * x, double * b) {  void Paso_Solver_solveAMG(Paso_Solver_AMG * amg, double * x, double * b) {
334       dim_t i,k;       dim_t i;
335       dim_t n_block=amg->n_block;       double time0=0;
336       double *r=MEMALLOC(amg->n,double);       double *r=NULL, *x0=NULL;
337       double *x0=MEMALLOC(amg->n,double);       bool_t verbose=0;
338         #ifdef MKL
339            Paso_SparseMatrix *temp=NULL;
340         #else
341            #ifdef UMFPACK
342               Paso_UMFPACK_Handler * ptr=NULL;
343            #endif
344         #endif
345         r=MEMALLOC(amg->n,double);
346         x0=MEMALLOC(amg->n,double);
347            
348       if (amg->level==0) {       if (amg->coarsest_level) {
349          /* x=invA_FF*b  */        
350          Paso_Solver_applyBlockDiagonalMatrix(n_block,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,x,b);        time0=Paso_timer();
351           #ifdef MKL
352            temp=Paso_SparseMatrix_alloc(MATRIX_FORMAT_BLK1 + MATRIX_FORMAT_OFFSET1, amg->A->pattern,1,1, FALSE);
353            #pragma omp parallel for private(i) schedule(static)
354            for (i=0;i<amg->A->len;++i) {
355                 temp->val[i]=amg->A->val[i];
356            }
357            Paso_MKL1(temp,x,b,0);
358            Paso_SparseMatrix_free(temp);
359           #else
360              #ifdef UMFPACK
361                 ptr=(Paso_UMFPACK_Handler *)(amg->solver);
362                 Paso_UMFPACK1(&ptr,amg->A,x,b,0);
363                 amg->solver=(void*) ptr;
364              #else
365                 Paso_Solver_solveJacobi(amg->GS,x,b);
366              #endif
367           #endif
368           time0=Paso_timer()-time0;
369           if (verbose) fprintf(stderr,"timing: DIRECT SOLVER: %e\n\n\n",time0);
370          
371       } else {       } else {
372           fprintf(stdout,"LEVEL %d \n",amg->level);          /* presmoothing */
373          /* presmoothing on (Shure, x, b, r) */           time0=Paso_timer();
374          /****************/           Paso_Solver_solveJacobi(amg->GS,x,b);
375           Paso_Solver_GS* GS=MEMALLOC(1,Paso_Solver_GS);           time0=Paso_timer()-time0;
376           GS=Paso_Solver_getGS(amg->A,-1);           if (verbose) fprintf(stderr,"timing: Presmooting: %e\n",time0);
377           Paso_Solver_solveGS(GS,x,b);          /* end of presmoothing */
378                    
379            
380             time0=Paso_timer();
381           #pragma omp parallel for private(i) schedule(static)           #pragma omp parallel for private(i) schedule(static)
382           for (i=0;i<amg->n;++i) r[i]=b[i];           for (i=0;i<amg->n;++i) r[i]=b[i];
383                    
384           /*r=b-Ax*/           /*r=b-Ax*/
385           Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);           Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
386         /****************/                
387          /* b->[b_F,b_C]     */          /* b->[b_F,b_C]     */
388          if (n_block==1) {          #pragma omp parallel for private(i) schedule(static)
389             #pragma omp parallel for private(i) schedule(static)          for (i=0;i<amg->n_F;++i) amg->b_F[i]=r[amg->rows_in_F[i]];
            for (i=0;i<amg->n_F;++i) amg->b_F[i]=r[amg->rows_in_F[i]];    /* was b istead of r */  
           #pragma omp parallel for private(i) schedule(static)  
            for (i=0;i<amg->n_C;++i) amg->b_C[i]=r[amg->rows_in_C[i]];  
         } else {  
            #pragma omp parallel for private(i,k) schedule(static)  
            for (i=0;i<amg->n_F;++i)  
                  for (k=0;k<n_block;k++) amg->b_F[amg->n_block*i+k]=r[n_block*amg->rows_in_F[i]+k];  
            #pragma omp parallel for private(i,k) schedule(static)  
            for (i=0;i<amg->n_C;++i)  
                  for (k=0;k<n_block;k++) amg->b_C[amg->n_block*i+k]=r[n_block*amg->rows_in_C[i]+k];  
         }  
   
         /****************/  
         /* Coursening part */  
390                    
391          /* r_F=A_FF^-1*r_F  */          #pragma omp parallel for private(i) schedule(static)
392          /*Paso_Solver_applyBlockDiagonalMatrix(n_block,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->b_F,amg->b_F);          for (i=0;i<amg->n_C;++i) amg->b_C[i]=r[amg->rows_in_C[i]];
         fprintf(stderr,"2\n");  
         */  
         /* r_C=r_C-A_CF*r_F */  
         /*Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_CF,amg->b_F,1.,amg->b_C);  
         fprintf(stderr,"3\n");  
         */  
        /****************/  
393    
           
394          /* x_F=invA_FF*b_F  */          /* x_F=invA_FF*b_F  */
395          Paso_Solver_applyBlockDiagonalMatrix(n_block,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);          Paso_Solver_applyBlockDiagonalMatrix(1,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);
396                    
397          /* b_C=b_C-A_CF*x_F */          /* b_C=b_C-A_CF*x_F */
398          Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_CF,amg->x_F,1.,amg->b_C);          Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_CF,amg->x_F,1.,amg->b_C);
399                    
400            time0=Paso_timer()-time0;
401            if (verbose) fprintf(stderr,"timing: Before next level: %e\n",time0);
402            
403          /* x_C=AMG(b_C)     */          /* x_C=AMG(b_C)     */
404          Paso_Solver_solveAMG(amg->AMG_of_Schur,amg->x_C,amg->b_C);          Paso_Solver_solveAMG(amg->AMG_of_Schur,amg->x_C,amg->b_C);
405                    
406            time0=Paso_timer();
407          /* b_F=b_F-A_FC*x_C */          /* b_F=b_F-A_FC*x_C */
408          Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_FC,amg->x_C,1.,amg->b_F);          Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_FC,amg->x_C,1.,amg->b_F);
409          /* x_F=invA_FF*b_F  */          /* x_F=invA_FF*b_F  */
410          Paso_Solver_applyBlockDiagonalMatrix(n_block,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);          Paso_Solver_applyBlockDiagonalMatrix(1,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);
411          /* x<-[x_F,x_C]     */          /* x<-[x_F,x_C]     */
412            
413          /* */          #pragma omp parallel for private(i) schedule(static)
414                    for (i=0;i<amg->n;++i) {
415                        if (amg->mask_C[i]>-1) {
416          if (n_block==1) {                   x[i]=amg->x_C[amg->mask_C[i]];
417             #pragma omp parallel for private(i) schedule(static)              } else {
418             for (i=0;i<amg->n;++i) {                   x[i]=amg->x_F[amg->mask_F[i]];
419                if (amg->mask_C[i]>-1) {              }
                   x[i]+=amg->x_C[amg->mask_C[i]];  
               } else {  
                   x[i]+=amg->x_F[amg->mask_F[i]];  
               }  
            }  
         } else {  
            #pragma omp parallel for private(i,k) schedule(static)  
            for (i=0;i<amg->n;++i) {  
                  if (amg->mask_C[i]>-1) {  
                      for (k=0;k<n_block;k++) x[n_block*i+k]+=amg->x_C[n_block*amg->mask_C[i]+k];  
                  } else {  
                      for (k=0;k<n_block;k++) x[n_block*i+k]+=amg->x_F[n_block*amg->mask_F[i]+k];  
                  }  
            }  
         }  
         /* all done */  
      /*post smoothing*/  
 /*     Paso_Solver_solveGS(GS,x0,b);  
      if (n_block==1) {  
            #pragma omp parallel for private(i) schedule(static)  
            for (i=0;i<amg->n;++i) x[i]+=x0[i];  
         } else {  
            #pragma omp parallel for private(i,k) schedule(static)  
            for (i=0;i<amg->n;++i) {  
                      for (k=0;k<n_block;k++) x[n_block*i+k]+=x0[n_block*i+k];  
            }  
420          }          }
421  */              
422       Paso_Solver_GS_free(GS);          time0=Paso_timer()-time0;
423            if (verbose) fprintf(stderr,"timing: After next level: %e\n",time0);
424    
425         /*postsmoothing*/
426         time0=Paso_timer();
427         #pragma omp parallel for private(i) schedule(static)
428         for (i=0;i<amg->n;++i) r[i]=b[i];
429        
430         /*r=b-Ax */
431         Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
432         Paso_Solver_solveJacobi(amg->GS,x0,r);
433        
434         #pragma omp parallel for private(i) schedule(static)
435         for (i=0;i<amg->n;++i) x[i]+=x0[i];
436        
437         time0=Paso_timer()-time0;
438         if (verbose) fprintf(stderr,"timing: Postsmoothing: %e\n",time0);
439    
440         /*end of postsmoothing*/
441        
442       }       }
      MEMFREE(x0);  
443       MEMFREE(r);       MEMFREE(r);
444         MEMFREE(x0);
445       return;       return;
446  }  }
   
 /*  
  * $Log$  
  *  
  */  

Legend:
Removed from v.1890  
changed lines
  Added in v.2551

  ViewVC Help
Powered by ViewVC 1.1.26