Actual source code: samgpetsc.c

  1: #include "petscksp.h"
 2:  #include src/mat/impls/aij/seq/aij.h
 3:  #include global.h
 4:  #include externc.h
 5:  #include samgfunc.h
 6:  #include petscfunc.h

  8: /*MC
  9:      PCSAMG -   SAMG + PETSc interface                                                
 10:                                                          
 11:     This interface allows e.g. to call samg as a shell preconditioner.    
 12:     samg is the new algebraic multigrid code by Klaus Stueben [1].         
 13:     Reference for SAMG                                                     
 14:     [1] K. St\"{u}ben,"Algebraic Multigrid: An Introduction for Positive  
 15:         Definite Problems with Applications", in [2], pp. 413--532,        
 16:         Academic Press, 2001.                                              
 17:     [2] U. Trottenberg, C. Oosterlee and A. Sch\"{u}ller, "Multigrid      
 18:         Methods", Academic Press, 2001.                                    
 19:     [1] is also available as                                               
 20:     [3] K. St\"{u}ben,"Algebraic Multigrid: An Introduction for Positive  
 21:         Definite Problems with Applications", Tech. Rep. 53, German       
 22:         National Research Center for Information Technology (GMD),        
 23:         Schloss Birlinhoven, D-53754 Sankt-Augustin, Germany, March 1999   
 24:     For more information on the SAMG-PETSc interface and examples of it's 
 25:     use, see                                                              
 26:     [4] D. Lahaye "Algebraic Multigrid for Time-Harmonic Magnetic Field    
 27:         Computations", PhD Thesis, KU Leuven, Belgium, December 2001.      
 28:         (www.cs.kuleuven.ac.be/~domenico)                                  
 29:                                                                           
 30:    Notes on PETSc part of this interface                                  
 31:    Information on how to set up shell preconditioners in PETSc can be     
 32:    in the PETSc documentation under preconditioners.                       

 34:    This preconditioner has not been completely organized to match the PETSc style,
 35:    see src/ksp/pc/impls/samgpetsc.c for the PC shell routines.

 37:     SAMG is a commercial product available from Klaus Stueben.
 38:  
 39:    Level: developer

 41:    Contributed by Domenico Lahaye (domenico.lahaye@cs.kuleuven.ac.be)   January 2001 
 42:                                                                           
 43: M*/

 45: /* ------------------------------------------------------------------- */
 48: /*.. SamgShellPCCreate - This routine creates a user-defined
 49:      preconditioner context.

 51:      Output Parameter:
 52:      shell - user-defined preconditioner context..*/

 54: PetscErrorCode SamgShellPCCreate(SamgShellPC **shell)
 55: {
 56:    SamgShellPC *newctx;

 59:    PetscNew(SamgShellPC,&newctx);
 60:    *shell = newctx;
 61:    return 0;
 62: }

 64: /* ------------------------------------------------------------------- */
 67: /*..SamgShellPCSetUp - This routine sets up a user-defined
 68:     ramg preconditioner context.  

 70:     Input Parameters:
 71:     shell - user-defined preconditioner context
 72:     pmat  - preconditioner matrix

 74:     Output Parameter:
 75:     shell - fully set up user-defined preconditioner context

 77:     This routine calls the setup phase of RAMG..*/
 78: 
 79: PetscErrorCode SamgShellPCSetUp(SamgShellPC *shell, Mat pmat)
 80: {
 82:    int  numnodes, numnonzero, nnz_count;
 83:    int              j, I, J, ncols_getrow, *cols_getrow, *diag;
 84:    PetscScalar      *vals_getrow;
 85:    MatInfo          info;
 86:    Mat_SeqAIJ       *aij = (Mat_SeqAIJ*)pmat->data;
 87:    /*..SAMG variables..*/
 88:    SAMG_PARAM       *samg_param;
 89:    double           *u_approx, *rhs, *Asky;
 90:    int              *ia, *ja;
 91:    /*..Primary SAMG parameters..*/
 92:    /*....System of equations....*/
 93:    int      matrix;
 94:    /*....Start solution and stopping criterium....*/
 95:    int      ifirst;
 96:    double   eps;
 97:    /*....Scalar and coupled system..*/
 98:    int      nsys=1, ndiu=1, ndip=1;
 99:    int      iscale[1], iu[1], ip[1];
100:    /*....Approach and smoother....*/
101:    int      nsolve;
102:    /*....Cycling process....*/
103:    int      ncyc;
104:    /*....Repeated calls....*/
105:    int      iswtch;
106:    /*....Initial dimensioning..*/
107:    double   a_cmplx, g_cmplx, p_cmplx, w_avrge;
108:    /*....Class 1 input parameters controlling input/output..*/
109:    int      idump, iout;
110:    double   chktol;
111:    /*....Output parameters....*/
112:    double   res_in, res_out;
113:    int      ncyc_done;
114:    /*..Numbers of levels created by SAMG..*/
115:    int      levels;
116:    /*..Auxilary integer to set secondary parameters..*/
117:    int      intin;

119:    /*..Get size and number of unknowns of preconditioner matrix..*/
120:    MatGetSize(pmat, &numnodes, &numnodes);
121:    MatGetInfo(pmat,MAT_LOCAL,&info);
122:    numnonzero = int(info.nz_used);

124:    /*..Allocate memory for RAMG variables..*/
125:    PetscMalloc(numnonzero   * sizeof(double),&Asky);
126:    PetscMalloc((numnodes+1) * sizeof(int),&ia);
127:    PetscMalloc(numnonzero   * sizeof(int),&ja);
128:    PetscMalloc(numnodes     * sizeof(double),&u_approx);
129:    PetscMalloc(numnodes     * sizeof(double),&rhs);

131:    /*..Store PETSc matrix in compressed skyline format required by SAMG..*/
132:    nnz_count = 0;
133:    MatMarkDiagonal_SeqAIJ(pmat);
134:    diag = aij->diag;

136:    for (I=0;I<numnodes;I++){
137:      ia[I]           = nnz_count;

139:      /*....put in diagonal entry first....*/
140:      ja[nnz_count]   = I;
141:      Asky[nnz_count] = aij->a[diag[I]];
142:      nnz_count++;

144:      /*....put in off-diagonals....*/
145:      ncols_getrow = aij->i[I+1] - aij->i[I];
146:      vals_getrow  = aij->a + aij->i[I];
147:      cols_getrow  = aij->j + aij->i[I];
148:      for (j=0;j<ncols_getrow;j++){
149:        J = cols_getrow[j];
150:        if (J != I) {
151:          Asky[nnz_count] = vals_getrow[j];
152:          ja[nnz_count]   = J;
153:          nnz_count++;
154:        }
155:      }
156:    }
157:    ia[numnodes] = nnz_count;

159:    /*..Allocate memory for SAMG parameters..*/
160:    PetscNew(SAMG_PARAM,&samg_param);

162:    /*..Set SAMG parameters..*/
163:    SamgGetParam(samg_param);

165:    /*..Set Primary parameters..*/
166:    matrix  = 12;
167:    ifirst  = (*samg_param).IFIRST;
168:    eps     = (*samg_param).EPS;
169:    nsolve  = (*samg_param).NSOLVE;
170:    ncyc    = (*samg_param).NCYC;
171:    iswtch  = (*samg_param).ISWTCH;
172:    a_cmplx = (*samg_param).A_CMPLX;
173:    g_cmplx = (*samg_param).G_CMPLX;
174:    p_cmplx = (*samg_param).P_CMPLX;
175:    w_avrge = (*samg_param).W_AVRGE;
176:    chktol  = (*samg_param).CHKTOL;
177:    idump   = (*samg_param).IDUMP;
178:    iout    = (*samg_param).IOUT;
179:    /*..Set secondary parameters..*/
180:    SAMG_set_levelx(&(samg_param->LEVELX));
181:    SAMG_set_nptmn(&(samg_param->NPTMN));
182:    SAMG_set_ecg(&(samg_param->ECG));
183:    SAMG_set_ewt(&(samg_param->EWT));
184:    SAMG_set_ncg(&(samg_param->NCG));
185:    SAMG_set_nwt(&(samg_param->NWT));
186:    SAMG_set_etr(&(samg_param->ETR));
187:    SAMG_set_ntr(&(samg_param->NTR));
188:    SAMG_set_nrd(&(samg_param->NRD));
189:    SAMG_set_nru(&(samg_param->NRU));
190:    SAMG_set_nrc(&(samg_param->NRC));
191:    intin = 6; SAMG_set_logio(&intin);
192:    intin = 1; SAMG_set_mode_debug(&intin);
193:    intin = 0; SAMG_set_iter_pre(&intin);
194:    intin = 0; SAMG_set_np_opt(&intin);
195:    intin = 0; SAMG_set_ncgrad_default(&intin);

197:    /*..Reset ncyc such that only setup is performed. This is done by setting 
198:      the last digit of ncyc (the number of cycles performed) equal to zero. 
199:      The first digits of ncyc (related to the solve phase) become 
200:      irrelevant.
201:    ..*/
202:    ncyc   = 1000;

204:    /*..Switch arrays ia and ja to Fortran conventions..*/
205:    for (j=0;j<=numnodes;j++)
206:        ia[j]++;
207:    for (j=0;j<numnonzero;j++)
208:        ja[j]++;

210:    PetscPrintf(MPI_COMM_WORLD,"\n\n");
211:    PetscPrintf(MPI_COMM_WORLD,"******************************************\n");
212:    PetscPrintf(MPI_COMM_WORLD,"*** Start Setup SAMG code (scal. mode) ***\n");
213:    PetscPrintf(MPI_COMM_WORLD,"******************************************\n");
214:    PetscPrintf(MPI_COMM_WORLD,"\n\n");

216:    /*..Call SAMG..*/
217:    SAMG(&numnodes, &numnonzero, &nsys,
218:         ia, ja, Asky, rhs, u_approx, iu, &ndiu, ip, &ndip, &matrix,
219:         iscale, &res_in, &res_out, &ncyc_done, &ierr, &nsolve,
220:         &ifirst, &eps, &ncyc, &iswtch, &a_cmplx, &g_cmplx,
221:          &p_cmplx, &w_avrge, &chktol, &idump, &iout);

223:    PetscPrintf(MPI_COMM_WORLD,"\n\n");
224:    PetscPrintf(MPI_COMM_WORLD,"******************************************\n");
225:    PetscPrintf(MPI_COMM_WORLD,"*** End Setup SAMG code (scal. mode)   ***\n");
226:    PetscPrintf(MPI_COMM_WORLD,"******************************************\n");
227:    PetscPrintf(MPI_COMM_WORLD,"\n\n");

229:    /*..Get number of levels created..*/
230:    SAMGPETSC_get_levels(&levels);

232:    /*..Store RAMG output in PETSc context..*/
233:    shell->A        = Asky;
234:    shell->IA       = ia;
235:    shell->JA       = ja;
236:    shell->PARAM    = samg_param;
237:    (*shell).LEVELS = levels;

239:    return 0;
240: }

242: /* ------------------------------------------------------------------- */
245: /*..SamgShellPCApply - This routine applies the AMG code as preconditioner 

247:     Input Parameters:
248:     ctx - user-defined context, as set by SamgShellSetApply()
249:     x - input vector

251:     Output Parameter:
252:     y - preconditioned vector

254:    Notes:
255:    Note that the PCSHELL preconditioner passes a void pointer as the
256:    first input argument.  This can be cast to be the whatever the user
257:    has set (via PCSetShellApply()) the application-defined context to be.

259:    ..*/
260: /*..To apply AMG as a preconditioner we set:                        */
261: /*        i) rhs-vector equal to the residual                       */
262: /*       ii) start solution equal to zero                           */
263: /*  Implementation notes:                                           */
264: /*  For the residual (vector r) we take the values from the vector  */
265: /*  using VecGetArray. No explicit memory allocation for            */
266: /*  vals_getarray is thus  needed as VecGetArray takes care of it.  */
267: /*  The allocated memory is freed again using a call to             */
268: /*  VecRestoreArray. The values in vals_getarray are then copy to   */
269: /*  rhs of the AMG code using memcpy.                               */

271: PetscErrorCode SamgShellPCApply(void *ctx, Vec r, Vec z)
272: {
274:    int  I, numnodes, numnonzero, *cols;
275:    SamgShellPC      *shell = (SamgShellPC *) ctx;
276:    double           *u_approx, *rhs, *Asky, *vals_getarray;
277:    int              *ia, *ja;
278:    SAMG_PARAM       *samg_param;
279:    /*..Primary SAMG parameters..*/
280:    /*....System of equations....*/
281:    int              matrix;
282:    /*....Start solution and stopping criterium....*/
283:    int              ifirst;
284:    double           eps;
285:    /*....Scalar and coupled system..*/
286:    int              nsys=1, ndiu=1, ndip=1;
287:    int              iscale[1], iu[1], ip[1];
288:    /*....Approach and smoother....*/
289:    int              nsolve;
290:    /*....Cycling process....*/
291:    int              ncyc;
292:    /*....Repeated calls....*/
293:    int              iswtch;
294:    /*....Initial dimensioning..*/
295:    double           a_cmplx, g_cmplx, p_cmplx, w_avrge;
296:    /*....Class 1 input parameters controlling input/output..*/
297:    int              idump, iout;
298:    double           chktol;
299:    /*....Output parameters....*/
300:    double           res_in, res_out;
301:    int              ncyc_done;

303:    /*..Get values from context..*/
304:    Asky       = shell->A;
305:    ia         = shell->IA;
306:    ja         = shell->JA;
307:    samg_param = shell->PARAM;

309:    /*..Get numnodes and numnonzeros..*/
310:    /*....numnodes can be determined as the size of the input vector r....*/
311:    VecGetSize(r,&numnodes);
312:    /*....numnonzero is determined from the pointer ia....*/
313:    /*....Remember that ia following Fortran conventions....*/
314:    numnonzero = ia[numnodes]-1;

316:    /*..Set the rhs of the call to ramg equal to the residual..*/
317:    VecGetArray(r,&vals_getarray);

319:    /*..Allocate memory for rhs and initial solution of call to samg..*/
320:    PetscMalloc(numnodes     * sizeof(double),&u_approx);
321:    PetscMalloc(numnodes     * sizeof(double),&rhs);

323:    /*..Set rhs of call to ramg..*/
324:    memcpy(rhs, vals_getarray, numnodes * sizeof(*rhs));
325: 
326:    /*..Set initial solution of call to ramg to zero..*/
327:    for (I=0;I<numnodes;I++){
328:        u_approx[I] = 0.;
329:    }

331:    /*..Set Primary parameters..*/
332:    matrix  = (*samg_param).MATRIX;
333:    ifirst  = (*samg_param).IFIRST;
334:    eps     = (*samg_param).EPS;
335:    nsolve  = (*samg_param).NSOLVE;
336:    ncyc    = (*samg_param).NCYC;
337:    iswtch  = (*samg_param).ISWTCH;
338:    a_cmplx = (*samg_param).A_CMPLX;
339:    g_cmplx = (*samg_param).G_CMPLX;
340:    p_cmplx = (*samg_param).P_CMPLX;
341:    w_avrge = (*samg_param).W_AVRGE;
342:    chktol  = (*samg_param).CHKTOL;
343:    idump   = (*samg_param).IDUMP;
344:    iout    = (*samg_param).IOUT;

346:    /*..Redefine iswtch to bypass setup..*/
347:    /*....First digit of iswtch = 2: bypass setup and do not release memory
348:          upon return....*/
349:    /*....Second digit of iswtch = 1: memory extension switch....*/
350:    /*....Third and fourth digit of iswtch: n_default. If n_default = 0, 
351:      the user has to set secondary parameters....*/
352:    iswtch = 210;

354:    /*..Call SAMG..*/
355:    SAMG(&numnodes, &numnonzero, &nsys,
356:         ia, ja, Asky, rhs, u_approx, iu, &ndiu, ip, &ndip, &matrix,
357:         iscale, &res_in, &res_out, &ncyc_done, &ierr, &nsolve,
358:         &ifirst, &eps, &ncyc, &iswtch, &a_cmplx, &g_cmplx,
359:         &p_cmplx, &w_avrge, &chktol, &idump, &iout);

361:    /*..Create auxilary integer array..*/
362:    PetscMalloc(numnodes * sizeof(double),&cols);

364:    for (I=0;I<numnodes;I++)
365:        cols[I] = I;

367:    /*..Store values computed by SAMG into the PETSc vector z..*/
368:    VecSetValues(z,numnodes,cols,u_approx,INSERT_VALUES);
369: 

371:    /*..Restore PETSc rhs vector..*/
372:    VecRestoreArray(r, &vals_getarray);

374:    PetscFree(cols);
375:    PetscFree(rhs);
376:    PetscFree(u_approx);

378:    return 0;
379: }

381: /* ------------------------------------------------------------------- */
384: /*..RamgShellPCDestroy - This routine destroys a user-defined
385:     preconditioner context.

387:     Input Parameter:
388:     shell - user-defined preconditioner context..*/

390: PetscErrorCode SamgShellPCDestroy(SamgShellPC *shell)
391: {
392:   /*..Free memory allocated by samg..*/
393:   SAMG_cleanup();
394:   /*..Free PCShell context..*/
395:   PetscFree(shell->A);
396:   PetscFree(shell->IA);
397:   PetscFree(shell->JA);
398:   PetscFree(shell->PARAM);
399:   PetscFree(shell);

401:   return 0;
402: }
403: /* ------------------------------------------------------------------- */
406: /*..SamgGetParam - Gets SAMG parameters specified at runtime 
407:     OUTPUT: The parameters set in the SAMG_PARAM context
408: ..*/

410: PetscErrorCode SamgGetParam(SAMG_PARAM *samg_param)
411: {

414:   /*..Set default SAMG paramets..*/
415:   /*....Class 0 SAMG parameters....*/
416:   (*samg_param).MATRIX    = 12;
417:   /*....Primary SAMG parameters....*/
418:   /*......If ifirst=0, the vector u, as passed to SAMG, is taken as first 
419:           approximation......*/
420:   (*samg_param).IFIRST    = 0;
421:   (*samg_param).EPS       = 1e-12;
422:   /*......nsolve =1 denotes scalar approach......*/
423:   (*samg_param).NSOLVE    = 1;
424:   /*......note: in the AMG-PETSc interface the number of cycles is required 
425:           to equal one to assure that in the PCApply routine AMG only performs 
426:           one cycle......*/
427:   (*samg_param).NCYC      = 1001;
428:   (*samg_param).ISWTCH    = 410;
429:   (*samg_param).A_CMPLX   = 2.5;
430:   (*samg_param).G_CMPLX   = 1.9;
431:   (*samg_param).P_CMPLX   = 1.9;
432:   (*samg_param).W_AVRGE   = 2.5;
433:   (*samg_param).CHKTOL    = 1e-8;
434:   (*samg_param).IDUMP     = 0;
435:   (*samg_param).IOUT      = 2;
436:   /*....Secundary SAMG parameters....*/
437:   (*samg_param).LEVELX    = 25;
438:   (*samg_param).NPTMN     = 100;
439:   (*samg_param).ECG       = 0.25;
440:   (*samg_param).EWT       = 0.5;
441:   (*samg_param).NCG       = 1000;
442:   (*samg_param).NWT       = 3000;
443:   (*samg_param).ETR       = 12.2;
444:   (*samg_param).NTR       = 2;
445:   (*samg_param).NRD       = 131;
446:   (*samg_param).NRC       = 0;
447:   (*samg_param).NRU       = -131;

449:   /*..Overwrite default values by values specified at runtime..*/
450:   /*....Primary SAMG parameters....*/
451:   PetscOptionsGetInt(PETSC_NULL,"-pc_samg_iswtch",&(*samg_param).ISWTCH,
452:                        PETSC_NULL);

454:   PetscOptionsGetInt(PETSC_NULL,"-pc_samg_ncyc",&(*samg_param).NCYC,
455:                        PETSC_NULL);

457:   PetscOptionsGetInt(PETSC_NULL,"-pc_samg_iout",&(*samg_param).IOUT,
458:                        PETSC_NULL);

460:  /*....Secundary SAMG parameters....*/
461:   PetscOptionsGetInt(PETSC_NULL,"-pc_samg_levelx",&(*samg_param).LEVELX,
462:                        PETSC_NULL);

464:   PetscOptionsGetInt(PETSC_NULL,"-pc_samg_nptmn",&(*samg_param).NPTMN,
465:                        PETSC_NULL);
466:   PetscOptionsGetInt(PETSC_NULL,"-pc_samg_nrd",&(*samg_param).NRD,
467:                        PETSC_NULL);

469:   PetscOptionsGetInt(PETSC_NULL,"-pc_samg_nrc",&(*samg_param).NRC,
470:                        PETSC_NULL);

472:   PetscOptionsGetInt(PETSC_NULL,"-pc_samg_nru",&(*samg_param).NRU,
473:                        PETSC_NULL);

475:   return 0;
476: }