# include # include # include # include # include # include "poisson_1d_multigrid.h" /******************************************************************************/ void ctof ( int nc, double uc[], int nf, double uf[] ) /******************************************************************************/ /* Purpose: ctof() transfers data from a coarse to a finer grid. Licensing: This code is distributed under the MIT license. Modified: 07 December 2011 Author: John Burkardt Reference: William Hager, Applied Numerical Linear Algebra, Prentice-Hall, 1988, ISBN13: 978-0130412942, LC: QA184.H33. Input: int NC, the number of coarse nodes. double UC[NC], the coarse correction data. int NF, the number of fine nodes. double UF[NF]: the current fine grid data.. Output: double UF[NF]: the updated fine grid data, using prolonged coarse correction data. */ { int ic; int iff; /* Add coarse node data to corresponding fine nodes. */ for ( ic = 0; ic < nc; ic++ ) { iff = 2 * ic; uf[iff] = uf[iff] + uc[ic]; } /* Add average of two coarse nodes to intermediate fine nodes. */ for ( ic = 0; ic < nc - 1; ic++ ) { iff = 2 * ic + 1; uf[iff] = uf[iff] + 0.5 * ( uc[ic] + uc[ic+1] ); } return; } /******************************************************************************/ void ftoc ( int nf, double uf[], double rf[], int nc, double uc[], double rc[] ) /******************************************************************************/ /* Purpose: ftoc() transfers data from a fine grid to a coarser grid. Licensing: This code is distributed under the MIT license. Modified: 07 December 2011 Author: John Burkardt Reference: William Hager, Applied Numerical Linear Algebra, Prentice-Hall, 1988, ISBN13: 978-0130412942, LC: QA184.H33. Input: int NF, the number of fine nodes. double UF[NF], the fine data. double RF[NF], the right hand side for the fine grid. int NC, the number of coarse nodes. Output: double UC[NC], the coarse grid data, set to zero. double RC[NC], the right hand side for the coarse grid. */ { int ic; int iff; for ( ic = 0; ic < nc; ic++ ) { uc[ic] = 0.0; } rc[0] = 0.0; for ( ic = 1; ic < nc - 1; ic++ ) { iff = 2 * ic; rc[ic] = 4.0 * ( rf[iff] + uf[iff-1] - 2.0 * uf[iff] + uf[iff+1] ); } rc[nc-1] = 0.0; return; } /******************************************************************************/ void gauss_seidel ( int n, double r[], double u[], double *dif_l1 ) /******************************************************************************/ /* Purpose: gauss_seidel() carries out one step of a Gauss-Seidel iteration. Licensing: This code is distributed under the MIT license. Modified: 07 December 2011 Author: John Burkardt Reference: William Hager, Applied Numerical Linear Algebra, Prentice-Hall, 1988, ISBN13: 978-0130412942, LC: QA184.H33. Input: int N, the number of unknowns. double R[N], the right hand side. double U[N], the estimated solution. Output: double U[N], the updated solution. double *DIF_L1, the L1 norm of the change in the solution estimates. */ { int i; double u_old; *dif_l1 = 0.0; for ( i = 1; i < n - 1; i++ ) { u_old = u[i]; u[i] = 0.5 * ( u[i-1] + u[i+1] + r[i] ); *dif_l1 = *dif_l1 + fabs ( u[i] - u_old ); } return; } /******************************************************************************/ int i4_log_2 ( int i ) /******************************************************************************/ /* Purpose: i4_log_2() returns the integer part of the logarithm base 2 of an I4. Example: I I4_LOG_10 ----- -------- 0 0 1 0 2 1 3 1 4 2 5 2 7 2 8 3 9 3 1000 9 1024 10 Discussion: I4_LOG_2 ( I ) + 1 is the number of binary digits in I. Licensing: This code is distributed under the MIT license. Modified: 23 October 2007 Author: John Burkardt Input: int I, the number whose logarithm base 2 is desired. Output: int I4_LOG_2, the integer part of the logarithm base 2 of the absolute value of X. */ { int i_abs; int two_pow; int value; if ( i == 0 ) { value = 0; } else { value = 0; two_pow = 2; i_abs = abs ( i ); while ( two_pow <= i_abs ) { value = value + 1; two_pow = two_pow * 2; } } return value; } /******************************************************************************/ void poisson_1d_multigrid ( int n, double a, double b, double ua, double ub, double force ( double x ), double exact ( double x ), int *it_num, double u[] ) /******************************************************************************/ /* Purpose: poisson_1d_multigrid() solves a 1D PDE using the multigrid method. Discussion: This routine solves a 1D boundary value problem of the form - U''(X) = F(X) for A < X < B, with boundary conditions U(A) = UA, U(B) = UB. The multigrid method is used. Licensing: This code is distributed under the MIT license. Modified: 26 July 2014 Author: Original FORTRAN77 version by William Hager. This version by John Burkardt. Reference: William Hager, Applied Numerical Linear Algebra, Prentice-Hall, 1988, ISBN13: 978-0130412942, LC: QA184.H33. Input: int N, the number of intervals. N must be a power of 2. double A, B, the left and right endpoints. double UA, UB, the left and right boundary values. double FORCE ( double x ): the function which evaluates the right hand side. double EXACT ( double x ): the function which evaluates the exact solution. Output: int *IT_NUM, the number of iterations. double U[N+1], the computed solution. */ { double d0; double d1; double h; int i; int it; int j; int k; int l; int ll; int m; int nl; double *r; double tol; double utol; double *uu; double *x; /* Determine if we have enough storage. */ k = i4_log_2 ( n ); if ( n != pow ( 2, k ) ) { printf ( "\n" ); printf ( "poisson_1d_multigrid(): Fatal error!\n" ); printf ( " N is not a power of 2.\n" ); exit ( 1 ); } nl = n + n + k - 2; /* Initialization. */ it = 4; *it_num = 0; tol = 0.0001; utol = 0.7; m = n; /* Set the nodes. */ x = r8vec_linspace_new ( n + 1, a, b ); /* Set the right hand side. */ r = ( double * ) malloc ( nl * sizeof ( double ) ); r[0] = ua; h = ( b - a ) / ( double ) ( n ); for ( i = 1; i < n; i++ ) { r[i] = h * h * force ( x[i] ); } r[n] = ub; uu = ( double * ) malloc ( nl * sizeof ( double ) ); for ( i = 0; i < nl; i++ ) { uu[i] = 0.0; } /* L points to the first entry of the solution. LL points to the penultimate entry. */ l = 0; ll = n - 1; /* Gauss-Seidel iteration */ d1 = 0.0; j = 0; for ( ; ; ) { d0 = d1; j = j + 1; gauss_seidel ( n + 1, r + l, uu + l, &d1 ); *it_num = *it_num + 1; /* Do at least 4 iterations at each level. */ if ( j < it ) { continue; } /* Enough iterations, satisfactory decrease, on finest grid, exit. */ else if ( d1 < tol && n == m ) { break; } /* Enough iterations, satisfactory convergence, go finer. */ else if ( d1 < tol ) { ctof ( n + 1, uu + l, n + n + 1, uu + (l-1-n-n) ); n = n + n; ll = l - 2; l = l - 1 - n; j = 0; } /* Enough iterations, slow convergence, 2 < N, go coarser. */ else if ( utol * d0 <= d1 && 2 < n ) { ftoc ( n + 1, uu + l, r + l, (n/2)+1, uu+(l+n+1), r+(l+n+1) ); n = n / 2; l = ll + 2; ll = ll + n + 1; j = 0; } } for ( i = 0; i < n + 1; i++ ) { u[i] = uu[i]; } /* Free memory. */ free ( r ); free ( uu ); free ( x ); return; } /******************************************************************************/ double *r8vec_linspace_new ( int n, double a, double b ) /******************************************************************************/ /* Purpose: r8vec_linspace_new() creates a vector of linearly spaced values. Discussion: An R8VEC is a vector of R8's. 4 points evenly spaced between 0 and 12 will yield 0, 4, 8, 12. In other words, the interval is divided into N-1 even subintervals, and the endpoints of intervals are used as the points. Licensing: This code is distributed under the MIT license. Modified: 29 March 2011 Author: John Burkardt Input: int N, the number of entries in the vector. double A, B, the first and last entries. Output: double R8VEC_LINSPACE_NEW[N], a vector of linearly spaced data. */ { int i; double *x; x = ( double * ) malloc ( n * sizeof ( double ) ); if ( n == 1 ) { x[0] = ( a + b ) / 2.0; } else { for ( i = 0; i < n; i++ ) { x[i] = ( ( double ) ( n - 1 - i ) * a + ( double ) ( i ) * b ) / ( double ) ( n - 1 ); } } return x; }