/[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 2625 by jfenwick, Fri Aug 21 06:30:25 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    
21  /* Author: artak@access.edu.au                                */  /* Author: artak@uq.edu.au                                */
22    
23  /**************************************************************/  /**************************************************************/
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);
45          Paso_SparseMatrix_free(in->A_CF);          Paso_SparseMatrix_free(in->A_CF);
46            Paso_SparseMatrix_free(in->A);
47          MEMFREE(in->rows_in_F);          MEMFREE(in->rows_in_F);
48          MEMFREE(in->rows_in_C);          MEMFREE(in->rows_in_C);
49          MEMFREE(in->mask_F);          MEMFREE(in->mask_F);
# Line 45  void Paso_Solver_AMG_free(Paso_Solver_AM Line 52  void Paso_Solver_AMG_free(Paso_Solver_AM
52          MEMFREE(in->b_F);          MEMFREE(in->b_F);
53          MEMFREE(in->x_C);          MEMFREE(in->x_C);
54          MEMFREE(in->b_C);          MEMFREE(in->b_C);
55            #ifdef MKL
56            #else
57               #ifdef UMFPACK
58               Paso_UMFPACK1_free((Paso_UMFPACK_Handler*)(in->solver));
59               #endif
60            #endif
61            in->solver=NULL;
62            MEMFREE(in->b_C);
63          MEMFREE(in);          MEMFREE(in);
64       }       }
65  }  }
# Line 68  to Line 83  to
83     then AMG is applied to S again until S becomes empty     then AMG is applied to S again until S becomes empty
84    
85  */  */
86  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) {
87    Paso_Solver_AMG* out=NULL;    Paso_Solver_AMG* out=NULL;
88      Paso_Pattern* temp1=NULL;
89      Paso_Pattern* temp2=NULL;
90      bool_t verbose=options->verbose;
91    dim_t n=A_p->numRows;    dim_t n=A_p->numRows;
92    dim_t n_block=A_p->row_block_size;    dim_t n_block=A_p->row_block_size;
93    index_t* mis_marker=NULL;      index_t* mis_marker=NULL;  
94    index_t* counter=NULL;    index_t* counter=NULL;
   double *rs=NULL;    
95    index_t iPtr,*index, *where_p, iPtr_s;    index_t iPtr,*index, *where_p, iPtr_s;
96    dim_t i,k,j,j0;    dim_t i,k,j;
97    Paso_SparseMatrix * schur=NULL;    Paso_SparseMatrix * schur=NULL;
98    Paso_SparseMatrix * schur_withFillIn=NULL;    Paso_SparseMatrix * schur_withFillIn=NULL;
99    schur_withFillIn=MEMALLOC(1,Paso_SparseMatrix);    double S=0;
100        char fname[6];
101        
102    double A11,A12,A13,A21,A22,A23,A31,A32,A33,D,time0,time1,time2;    /*Make sure we have block sizes 1*/
103        if (A_p->col_block_size>1) {
104         Paso_setError(TYPE_ERROR,"Paso_Solver_getAMG: AMG requires column block size 1.");
105         return NULL;
106      }
107      if (A_p->row_block_size>1) {
108         Paso_setError(TYPE_ERROR,"Paso_Solver_getAMG: AMG requires row block size 1.");
109         return NULL;
110      }
111      out=MEMALLOC(1,Paso_Solver_AMG);
112    /* identify independend set of rows/columns */    /* identify independend set of rows/columns */
113    mis_marker=TMPMEMALLOC(n,index_t);    mis_marker=TMPMEMALLOC(n,index_t);
114    counter=TMPMEMALLOC(n,index_t);    counter=TMPMEMALLOC(n,index_t);
115    rs=TMPMEMALLOC(n,double);    if ( !( Paso_checkPtr(mis_marker) || Paso_checkPtr(counter) || Paso_checkPtr(out)) ) {
116    out=MEMALLOC(1,Paso_Solver_AMG);       out->AMG_of_Schur=NULL;
117    out->AMG_of_Schur=NULL;       out->inv_A_FF=NULL;
118    out->inv_A_FF=NULL;       out->A_FF_pivot=NULL;
119    out->A_FF_pivot=NULL;       out->A_FC=NULL;
120    out->A_FC=NULL;       out->A_CF=NULL;
121    out->A_CF=NULL;       out->rows_in_F=NULL;
122    out->rows_in_F=NULL;       out->rows_in_C=NULL;
123    out->rows_in_C=NULL;       out->mask_F=NULL;
124    out->mask_F=NULL;       out->mask_C=NULL;
125    out->mask_C=NULL;       out->x_F=NULL;
126    out->x_F=NULL;       out->b_F=NULL;
127    out->b_F=NULL;       out->x_C=NULL;
128    out->x_C=NULL;       out->b_C=NULL;
129    out->b_C=NULL;       out->GS=NULL;
130    out->A=Paso_SparseMatrix_getReference(A_p);       out->A=Paso_SparseMatrix_getReference(A_p);
131    out->level=level;       out->GS=NULL;
132         out->solver=NULL;
133         /*out->GS=Paso_Solver_getGS(A_p,verbose);*/
134         out->level=level;
135        
136           sprintf(fname,"A%d.mat",level);
137        
138   /* fprintf(stderr,"START OF MATRIX \n\n");         Paso_SparseMatrix_saveMM(A_p,fname);
139    for (i = 0; i < A_p->numRows; ++i) {    
140       for (iPtr=A_p->pattern->ptr[i];iPtr<A_p->pattern->ptr[i + 1]; ++iPtr) {       if (level==0 || n<=options->min_coarse_matrix_size) {
141         j=A_p->pattern->index[iPtr];           out->coarsest_level=TRUE;
142         fprintf(stderr,"A[%d,%d]=%.2f ",i,j,A_p->val[iPtr]);           #ifdef MKL
143       }           #else
144       fprintf(stderr,"\n");              #ifdef UMFPACK
145     }              #else
146     fprintf(stderr,"END OF MATRIX \n\n");                  out->GS=Paso_Solver_getJacobi(A_p);
147   */              #endif
148    if ( !(Paso_checkPtr(mis_marker) || Paso_checkPtr(out) || Paso_checkPtr(counter) ) ) {           #endif
149       /* identify independend set of rows/columns */       } else {
150       time0=Paso_timer();           out->coarsest_level=FALSE;
151       #pragma omp parallel for private(i) schedule(static)           out->GS=Paso_Solver_getJacobi(A_p);
152       for (i=0;i<n;++i) mis_marker[i]=-1;  
153       Paso_Pattern_RS(A_p,mis_marker,0.25);           /* identify independend set of rows/columns */
154       /*           #pragma omp parallel for private(i) schedule(static)
155       for (i=0;i<n;++i) fprintf(stderr," i=%d mis[i]=%d \n",i,mis_marker[i]);           for (i=0;i<n;++i) mis_marker[i]=-1;
156       */  
157       time2=Paso_timer()-time0;           if (options->coarsening_method == PASO_YAIR_SHAPIRA_COARSENING) {
158       if (Paso_noError()) {                Paso_Pattern_coup(A_p,mis_marker,options->coarsening_threshold);
159          #pragma omp parallel for private(i) schedule(static)           }
160          for (i = 0; i < n; ++i) counter[i]=mis_marker[i];           else if (options->coarsening_method == PASO_RUGE_STUEBEN_COARSENING) {
161          out->n=n;                Paso_Pattern_RS(A_p,mis_marker,options->coarsening_threshold);
162          out->n_block=n_block;           }
163          out->n_F=Paso_Util_cumsum(n,counter);           else if (options->coarsening_method == PASO_AGGREGATION_COARSENING) {
164          out->mask_F=MEMALLOC(n,index_t);               Paso_Pattern_Aggregiation(A_p,mis_marker,options->coarsening_threshold);
165          out->rows_in_F=MEMALLOC(out->n_F,index_t);          }
166          out->inv_A_FF=MEMALLOC(n_block*n_block*out->n_F,double);          else {
167          out->A_FF_pivot=NULL; /* later use for block size>3 */             /*Default coarseneing*/
168          if (! (Paso_checkPtr(out->mask_F) || Paso_checkPtr(out->inv_A_FF) || Paso_checkPtr(out->rows_in_F) ) ) {              Paso_Pattern_RS(A_p,mis_marker,options->coarsening_threshold);
169             #pragma omp parallel  
170             {               /*Paso_Pattern_Aggregiation(A_p,mis_marker,options->coarsening_threshold);*/
171            }
172        
173            if (Paso_noError()) {
174               #pragma omp parallel for private(i) schedule(static)
175               for (i = 0; i < n; ++i) counter[i]=mis_marker[i];
176    
177               out->n=n;
178               out->n_block=n_block;
179               out->n_F=Paso_Util_cumsum(n,counter);
180               out->mask_F=MEMALLOC(n,index_t);
181               out->rows_in_F=MEMALLOC(out->n_F,index_t);
182               out->inv_A_FF=MEMALLOC(n_block*n_block*out->n_F,double);
183               out->A_FF_pivot=NULL; /* later use for block size>3 */
184               if (! (Paso_checkPtr(out->mask_F) || Paso_checkPtr(out->inv_A_FF) || Paso_checkPtr(out->rows_in_F) ) ) {
185                /* creates an index for F from mask */                /* creates an index for F from mask */
186                #pragma omp for private(i) schedule(static)                #pragma omp parallel for private(i) schedule(static)
187                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;
188                #pragma omp for private(i) schedule(static)                #pragma omp parallel for private(i) schedule(static)
189                for (i = 0; i < n; ++i) {                for (i = 0; i < n; ++i) {
190                   if  (mis_marker[i]) {                   if  (mis_marker[i]) {
191                          out->rows_in_F[counter[i]]=i;                          out->rows_in_F[counter[i]]=i;
# Line 151  Paso_Solver_AMG* Paso_Solver_getAMG(Paso Line 194  Paso_Solver_AMG* Paso_Solver_getAMG(Paso
194                          out->mask_F[i]=-1;                          out->mask_F[i]=-1;
195                   }                   }
196                }                }
197                /* Compute row-sum for getting rs(A_FF)*/                
198                #pragma omp for private(i,iPtr) schedule(static)                /* Compute row-sum for getting rs(A_FF)^-1*/
199                  #pragma omp parallel for private(i,iPtr,j,S) schedule(static)
200                for (i = 0; i < out->n_F; ++i) {                for (i = 0; i < out->n_F; ++i) {
201                  rs[i]=0;                  S=0;
202    /*printf("%d : %d -> ",i, out->rows_in_F[i]);*/
203                  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) {
204                   rs[i]+=A_p->val[iPtr];                   j=A_p->pattern->index[iPtr];
205    /*if (j==out->rows_in_F[i]) printf("diagonal %e",A_p->val[iPtr]);*/
206                     if (mis_marker[j] && j==out->rows_in_F[i])
207                         S=A_p->val[iPtr];
208                  }                  }
209    /*printf("-> %e\n",S);*/
210                    out->inv_A_FF[i]=1./S;
211                }                }
212                           }
213                #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 */  
214    
215             if( Paso_noError()) {          if ( Paso_noError()) {
216                /* 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 */
217                out->n_C=n-out->n_F;                out->n_C=n-out->n_F;
218                 /*if (out->n_C>11) {*/              
219                 if (level>0) {               /*if (out->n_F>500) */
220                     out->rows_in_C=MEMALLOC(out->n_C,index_t);                out->rows_in_C=MEMALLOC(out->n_C,index_t);
221                     out->mask_C=MEMALLOC(n,index_t);                out->mask_C=MEMALLOC(n,index_t);
222                     if (! (Paso_checkPtr(out->mask_C) || Paso_checkPtr(out->rows_in_C) ) ) {                if (! (Paso_checkPtr(out->mask_C) || Paso_checkPtr(out->rows_in_C) ) ) {
223                         /* creates an index for C from mask */                     /* creates an index for C from mask */
224                         #pragma omp parallel for private(i) schedule(static)                     #pragma omp parallel for private(i) schedule(static)
225                         for (i = 0; i < n; ++i) counter[i]=! mis_marker[i];                     for (i = 0; i < n; ++i) counter[i]=! mis_marker[i];
226                         Paso_Util_cumsum(n,counter);                     Paso_Util_cumsum(n,counter);
227                         #pragma omp parallel                     #pragma omp parallel for private(i) schedule(static)
228                         {                     for (i = 0; i < out->n_C; ++i) out->rows_in_C[i]=-1;
229                            #pragma omp for private(i) schedule(static)                     #pragma omp parallel for private(i) schedule(static)
230                            for (i = 0; i < out->n_C; ++i) out->rows_in_C[i]=-1;                     for (i = 0; i < n; ++i) {
231                            #pragma omp for private(i) schedule(static)                              if  (! mis_marker[i]) {
                           for (i = 0; i < n; ++i) {  
                              if  (! mis_marker[i]) {  
232                                  out->rows_in_C[counter[i]]=i;                                  out->rows_in_C[counter[i]]=i;
233                                  out->mask_C[i]=counter[i];                                  out->mask_C[i]=counter[i];
234                               } else {                               } else {
235                                  out->mask_C[i]=-1;                                  out->mask_C[i]=-1;
236                               }                               }
237                            }                     }
238                        } /* end parallel region */                }
239                        /* get A_CF block: */          }
240                        out->A_CF=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_F,out->rows_in_C,out->mask_F);          if ( Paso_noError()) {
241                        if (Paso_noError()) {                  /* get A_CF block: */
242                           /* 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);
243                           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: */
244                           /* 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);
245                           if (Paso_noError()) {                  /* get A_CC block: */
246                              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);
247                                
248                              /*find the pattern of the schur complement with fill in*/          }
249                              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()) {
250                                               /*find the pattern of the schur complement with fill in*/
251                              /* copy values over*/                temp1=Paso_Pattern_multiply(PATTERN_FORMAT_DEFAULT,out->A_CF->pattern,out->A_FC->pattern);
252                              #pragma omp for private(i,iPtr,iPtr_s,j,j0) schedule(static)                temp2=Paso_Pattern_binop(PATTERN_FORMAT_DEFAULT, schur->pattern, temp1);
253                              for (i = 0; i < schur_withFillIn->numRows; ++i) {                schur_withFillIn=Paso_SparseMatrix_alloc(A_p->type,temp2,1,1, TRUE);
254                                for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {                Paso_Pattern_free(temp1);
255                                  j=schur_withFillIn->pattern->index[iPtr];                Paso_Pattern_free(temp2);
256                                  schur_withFillIn->val[iPtr]=0.;          }
257                                  for (iPtr_s=schur->pattern->ptr[i];iPtr_s<schur->pattern->ptr[i + 1]; ++iPtr_s){          if ( Paso_noError()) {
258                                      j0=schur->pattern->index[iPtr_s];                /* copy values over*/
259                                      if (j==j0) {                #pragma omp parallel for private(i,iPtr,j,iPtr_s,index,where_p) schedule(static)
260                                        schur_withFillIn->val[iPtr]=schur->val[iPtr_s];                for (i = 0; i < schur_withFillIn->numRows; ++i) {
261                                        break;                  for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {
262                                      }                     j=schur_withFillIn->pattern->index[iPtr];
263                                  }                     iPtr_s=schur->pattern->ptr[i];
264                                }                     schur_withFillIn->val[iPtr]=0.;
265                              }                     index=&(schur->pattern->index[iPtr_s]);
266                                                   where_p=(index_t*)bsearch(&j,
267                            /*  for (i = 0; i < schur_withFillIn->numRows; ++i) {                                          index,
268                                for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {                                          schur->pattern->ptr[i + 1]-schur->pattern->ptr[i],
269                                  j=schur_withFillIn->pattern->index[iPtr];                                          sizeof(index_t),
270                                  fprintf(stderr,"A_CC[%d,%d]=%.2f ",i,j,schur_withFillIn->val[iPtr]);                                          Paso_comparIndex);
271                                }                     if (where_p!=NULL) {
272                                fprintf(stderr,"\n");                            schur_withFillIn->val[iPtr]=schur->val[iPtr_s+(index_t)(where_p-index)];
273                              }*/                     }
                             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 */  
                               }  
                             }  
                          }  
                      }  
274                   }                   }
275                }                }
276             }                Paso_Solver_updateIncompleteSchurComplement(schur_withFillIn,out->A_CF,out->inv_A_FF,out->A_FF_pivot,out->A_FC);
277                  out->AMG_of_Schur=Paso_Solver_getAMG(schur_withFillIn,level-1,options);
278            }
279            /* allocate work arrays for AMG application */
280            if (Paso_noError()) {
281                       out->x_F=MEMALLOC(n_block*out->n_F,double);
282                       out->b_F=MEMALLOC(n_block*out->n_F,double);
283                       out->x_C=MEMALLOC(n_block*out->n_C,double);
284                       out->b_C=MEMALLOC(n_block*out->n_C,double);
285    
286                       if (! (Paso_checkPtr(out->x_F) || Paso_checkPtr(out->b_F) || Paso_checkPtr(out->x_C) || Paso_checkPtr(out->b_C) ) ) {
287                           #pragma omp parallel for private(i,k) schedule(static)
288                           for (i = 0; i < out->n_F; ++i) {
289                                for (k=0; k<n_block;++k) {
290                                       out->x_F[i*n_block+k]=0.;
291                                       out->b_F[i*n_block+k]=0.;
292                                }
293                            }
294                            #pragma omp parallel for private(i,k) schedule(static)
295                            for (i = 0; i < out->n_C; ++i) {
296                               for (k=0; k<n_block;++k) {
297                                   out->x_C[i*n_block+k]=0.;
298                                   out->b_C[i*n_block+k]=0.;
299                               }
300                            }
301                       }
302          }          }
303       }       }
304    }    }
305    TMPMEMFREE(mis_marker);    TMPMEMFREE(mis_marker);
306    TMPMEMFREE(counter);    TMPMEMFREE(counter);
307    TMPMEMFREE(rs);    Paso_SparseMatrix_free(schur);
308      Paso_SparseMatrix_free(schur_withFillIn);
309    if (Paso_noError()) {    if (Paso_noError()) {
310        if (verbose) {        if (verbose && level>0 && !out->coarsest_level) {
311           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);  
          }  
312       }       }
313       return out;       return out;
314    } else  {    } else  {
# Line 320  Paso_Solver_AMG* Paso_Solver_getAMG(Paso Line 323  Paso_Solver_AMG* Paso_Solver_getAMG(Paso
323    
324       in fact it solves       in fact it solves
325    
326    [      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]
327    [ 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]
328    
329   in the form   in the form
# Line 339  Paso_Solver_AMG* Paso_Solver_getAMG(Paso Line 342  Paso_Solver_AMG* Paso_Solver_getAMG(Paso
342  */  */
343    
344  void Paso_Solver_solveAMG(Paso_Solver_AMG * amg, double * x, double * b) {  void Paso_Solver_solveAMG(Paso_Solver_AMG * amg, double * x, double * b) {
345       dim_t i,k;       dim_t i;
346       dim_t n_block=amg->n_block;       double time0=0;
347       double *r=MEMALLOC(amg->n,double);       double *r=NULL, *x0=NULL;
348       double *x0=MEMALLOC(amg->n,double);       bool_t verbose=0;
349         #ifdef MKL
350            Paso_SparseMatrix *temp=NULL;
351         #else
352            #ifdef UMFPACK
353               Paso_UMFPACK_Handler * ptr=NULL;
354            #endif
355         #endif
356         r=MEMALLOC(amg->n,double);
357         x0=MEMALLOC(amg->n,double);
358            
359       if (amg->level==0) {       if (amg->coarsest_level) {
360          /* x=invA_FF*b  */        
361          Paso_Solver_applyBlockDiagonalMatrix(n_block,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,x,b);        time0=Paso_timer();
362           #ifdef MKL
363            temp=Paso_SparseMatrix_alloc(MATRIX_FORMAT_BLK1 + MATRIX_FORMAT_OFFSET1, amg->A->pattern,1,1, FALSE);
364            #pragma omp parallel for private(i) schedule(static)
365            for (i=0;i<amg->A->len;++i) {
366                 temp->val[i]=amg->A->val[i];
367            }
368            Paso_MKL1(temp,x,b,0);
369            Paso_SparseMatrix_free(temp);
370           #else
371              #ifdef UMFPACK
372                 ptr=(Paso_UMFPACK_Handler *)(amg->solver);
373                 Paso_UMFPACK1(&ptr,amg->A,x,b,0);
374                 amg->solver=(void*) ptr;
375              #else
376                 Paso_Solver_solveJacobi(amg->GS,x,b);
377              #endif
378           #endif
379           time0=Paso_timer()-time0;
380           if (verbose) fprintf(stderr,"timing: DIRECT SOLVER: %e\n\n\n",time0);
381          
382       } else {       } else {
383           fprintf(stdout,"LEVEL %d \n",amg->level);          /* presmoothing */
384          /* presmoothing on (Shure, x, b, r) */           time0=Paso_timer();
385          /****************/           Paso_Solver_solveJacobi(amg->GS,x,b);
386           Paso_Solver_GS* GS=MEMALLOC(1,Paso_Solver_GS);           time0=Paso_timer()-time0;
387           GS=Paso_Solver_getGS(amg->A,-1);           if (verbose) fprintf(stderr,"timing: Presmooting: %e\n",time0);
388           Paso_Solver_solveGS(GS,x,b);          /* end of presmoothing */
389            
390                    
391             time0=Paso_timer();
392           #pragma omp parallel for private(i) schedule(static)           #pragma omp parallel for private(i) schedule(static)
393           for (i=0;i<amg->n;++i) r[i]=b[i];           for (i=0;i<amg->n;++i) r[i]=b[i];
394                    
395           /*r=b-Ax*/           /*r=b-Ax*/
396           Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);           Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
397         /****************/                
398          /* b->[b_F,b_C]     */          /* b->[b_F,b_C]     */
399          if (n_block==1) {          #pragma omp parallel for private(i) schedule(static)
400             #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 */  
401                    
402          /* r_F=A_FF^-1*r_F  */          #pragma omp parallel for private(i) schedule(static)
403          /*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");  
         */  
        /****************/  
404    
           
405          /* x_F=invA_FF*b_F  */          /* x_F=invA_FF*b_F  */
406          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);
407                    
408          /* b_C=b_C-A_CF*x_F */          /* b_C=b_C-A_CF*x_F */
409          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);
410                    
411            time0=Paso_timer()-time0;
412            if (verbose) fprintf(stderr,"timing: Before next level: %e\n",time0);
413            
414          /* x_C=AMG(b_C)     */          /* x_C=AMG(b_C)     */
415          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);
416                    
417            time0=Paso_timer();
418          /* b_F=b_F-A_FC*x_C */          /* b_F=b_F-A_FC*x_C */
419          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);
420          /* x_F=invA_FF*b_F  */          /* x_F=invA_FF*b_F  */
421          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);
422          /* x<-[x_F,x_C]     */          /* x<-[x_F,x_C]     */
423            
424          /* */          #pragma omp parallel for private(i) schedule(static)
425                    for (i=0;i<amg->n;++i) {
426                        if (amg->mask_C[i]>-1) {
427          if (n_block==1) {                   x[i]=amg->x_C[amg->mask_C[i]];
428             #pragma omp parallel for private(i) schedule(static)              } else {
429             for (i=0;i<amg->n;++i) {                   x[i]=amg->x_F[amg->mask_F[i]];
430                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];  
            }  
431          }          }
432  */              
433       Paso_Solver_GS_free(GS);          time0=Paso_timer()-time0;
434            if (verbose) fprintf(stderr,"timing: After next level: %e\n",time0);
435    
436         /*postsmoothing*/
437         time0=Paso_timer();
438         #pragma omp parallel for private(i) schedule(static)
439         for (i=0;i<amg->n;++i) r[i]=b[i];
440        
441         /*r=b-Ax */
442         Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
443         Paso_Solver_solveJacobi(amg->GS,x0,r);
444        
445         #pragma omp parallel for private(i) schedule(static)
446         for (i=0;i<amg->n;++i) x[i]+=x0[i];
447        
448         time0=Paso_timer()-time0;
449         if (verbose) fprintf(stderr,"timing: Postsmoothing: %e\n",time0);
450    
451         /*end of postsmoothing*/
452        
453       }       }
      MEMFREE(x0);  
454       MEMFREE(r);       MEMFREE(r);
455         MEMFREE(x0);
456       return;       return;
457  }  }
   
 /*  
  * $Log$  
  *  
  */  

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

  ViewVC Help
Powered by ViewVC 1.1.26