xref: /petsc/src/snes/tutorials/ex31.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
1fcf85c8cSAdelekeBankole 
2fcf85c8cSAdelekeBankole static char help[] = "A Chebyshev spectral method for the compressible Blasius boundary layer equations.\n\n";
3fcf85c8cSAdelekeBankole 
4fcf85c8cSAdelekeBankole /*
5fcf85c8cSAdelekeBankole    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
6fcf85c8cSAdelekeBankole    file automatically includes:
7fcf85c8cSAdelekeBankole      petscsys.h       - base PETSc routines   petscvec.h - vectors
8fcf85c8cSAdelekeBankole      petscmat.h - matrices
9fcf85c8cSAdelekeBankole      petscis.h     - index sets            petscksp.h - Krylov subspace methods
10fcf85c8cSAdelekeBankole      petscviewer.h - viewers               petscpc.h  - preconditioners
11fcf85c8cSAdelekeBankole      petscksp.h   - linear solvers
12fcf85c8cSAdelekeBankole    Include "petscdt.h" so that we can have support for use of Quadrature formulas
13fcf85c8cSAdelekeBankole */
14fcf85c8cSAdelekeBankole /*F
15fcf85c8cSAdelekeBankole This examples solves the compressible Blasius boundary layer equations
16fcf85c8cSAdelekeBankole 2(\rho\muf'')' + ff'' = 0
17fcf85c8cSAdelekeBankole (\rho\muh')' + Prfh' + Pr(\gamma-1)Ma^{2}\rho\muf''^{2} = 0
18fcf85c8cSAdelekeBankole following Howarth-Dorodnitsyn transformation with boundary conditions
19fcf85c8cSAdelekeBankole f(0) = f'(0) = 0, f'(\infty) = 1, h(\infty) = 1, h = \theta(0). Where \theta = T/T_{\infty}
20fcf85c8cSAdelekeBankole Note: density (\rho) and viscosity (\mu) are treated as constants in this example
21fcf85c8cSAdelekeBankole F*/
22fcf85c8cSAdelekeBankole #include <petscsnes.h>
23fcf85c8cSAdelekeBankole #include <petscdt.h>
24fcf85c8cSAdelekeBankole 
25fcf85c8cSAdelekeBankole /*
26fcf85c8cSAdelekeBankole    User-defined routines
27fcf85c8cSAdelekeBankole */
28fcf85c8cSAdelekeBankole 
29fcf85c8cSAdelekeBankole extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *);
30fcf85c8cSAdelekeBankole 
31fcf85c8cSAdelekeBankole typedef struct {
32fcf85c8cSAdelekeBankole   PetscReal  Ma, Pr, h_0;
33fcf85c8cSAdelekeBankole   PetscInt   N;
34fcf85c8cSAdelekeBankole   PetscReal  dx_deta;
35fcf85c8cSAdelekeBankole   PetscReal *x;
36fcf85c8cSAdelekeBankole   PetscReal  gamma;
37fcf85c8cSAdelekeBankole } Blasius;
38fcf85c8cSAdelekeBankole 
399371c9d4SSatish Balay int main(int argc, char **argv) {
40fcf85c8cSAdelekeBankole   SNES        snes; /* nonlinear solver context */
41fcf85c8cSAdelekeBankole   Vec         x, r; /* solution, residual vectors */
42fcf85c8cSAdelekeBankole   PetscMPIInt size;
43fcf85c8cSAdelekeBankole   Blasius    *blasius;
44fcf85c8cSAdelekeBankole   PetscReal   L, *weight; /* L is size of the domain */
45fcf85c8cSAdelekeBankole 
46327415f7SBarry Smith   PetscFunctionBeginUser;
47fcf85c8cSAdelekeBankole   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
48fcf85c8cSAdelekeBankole   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
49fcf85c8cSAdelekeBankole   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Example is only for sequential runs");
50fcf85c8cSAdelekeBankole 
51fcf85c8cSAdelekeBankole   // Read command-line arguments
52fcf85c8cSAdelekeBankole   PetscCall(PetscCalloc1(1, &blasius));
53fcf85c8cSAdelekeBankole   blasius->Ma    = 2;   /* Mach number */
54fcf85c8cSAdelekeBankole   blasius->Pr    = 0.7; /* Prandtl number */
55fcf85c8cSAdelekeBankole   blasius->h_0   = 2.;  /* relative temperature at the wall */
56fcf85c8cSAdelekeBankole   blasius->N     = 10;  /* Number of Chebyshev terms */
57fcf85c8cSAdelekeBankole   blasius->gamma = 1.4; /* specific heat ratio */
58fcf85c8cSAdelekeBankole   L              = 5;
59fcf85c8cSAdelekeBankole   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Compressible Blasius boundary layer equations", "");
60fcf85c8cSAdelekeBankole   PetscCall(PetscOptionsReal("-mach", "Mach number at freestream", "", blasius->Ma, &blasius->Ma, NULL));
61fcf85c8cSAdelekeBankole   PetscCall(PetscOptionsReal("-prandtl", "Prandtl number", "", blasius->Pr, &blasius->Pr, NULL));
62fcf85c8cSAdelekeBankole   PetscCall(PetscOptionsReal("-h_0", "Relative enthalpy at wall", "", blasius->h_0, &blasius->h_0, NULL));
63fcf85c8cSAdelekeBankole   PetscCall(PetscOptionsReal("-gamma", "Ratio of specific heats", "", blasius->gamma, &blasius->gamma, NULL));
64fcf85c8cSAdelekeBankole   PetscCall(PetscOptionsInt("-N", "Number of Chebyshev terms for f", "", blasius->N, &blasius->N, NULL));
65fcf85c8cSAdelekeBankole   PetscCall(PetscOptionsReal("-L", "Extent of the domain", "", L, &L, NULL));
66fcf85c8cSAdelekeBankole   PetscOptionsEnd();
67fcf85c8cSAdelekeBankole   blasius->dx_deta = 2 / L; /* this helps to map [-1,1] to [0,L] */
68fcf85c8cSAdelekeBankole   PetscCall(PetscMalloc2(blasius->N - 3, &blasius->x, blasius->N - 3, &weight));
69fcf85c8cSAdelekeBankole   PetscCall(PetscDTGaussQuadrature(blasius->N - 3, -1., 1., blasius->x, weight));
70fcf85c8cSAdelekeBankole 
71fcf85c8cSAdelekeBankole   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72fcf85c8cSAdelekeBankole      Create nonlinear solver context
73fcf85c8cSAdelekeBankole      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74fcf85c8cSAdelekeBankole   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
75fcf85c8cSAdelekeBankole 
76fcf85c8cSAdelekeBankole   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77fcf85c8cSAdelekeBankole      Create matrix and vector data structures; set corresponding routines
78fcf85c8cSAdelekeBankole      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79fcf85c8cSAdelekeBankole   /*
80fcf85c8cSAdelekeBankole      Create vectors for solution and nonlinear function
81fcf85c8cSAdelekeBankole   */
82fcf85c8cSAdelekeBankole   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
83fcf85c8cSAdelekeBankole   PetscCall(VecSetSizes(x, PETSC_DECIDE, 2 * blasius->N - 1));
84fcf85c8cSAdelekeBankole   PetscCall(VecSetFromOptions(x));
85fcf85c8cSAdelekeBankole   PetscCall(VecDuplicate(x, &r));
86fcf85c8cSAdelekeBankole 
87fcf85c8cSAdelekeBankole   /*
88fcf85c8cSAdelekeBankole       Set function evaluation routine and vector.
89fcf85c8cSAdelekeBankole    */
90fcf85c8cSAdelekeBankole   PetscCall(SNESSetFunction(snes, r, FormFunction, blasius));
91fcf85c8cSAdelekeBankole   {
92fcf85c8cSAdelekeBankole     KSP ksp;
93fcf85c8cSAdelekeBankole     PC  pc;
94fcf85c8cSAdelekeBankole     SNESGetKSP(snes, &ksp);
95fcf85c8cSAdelekeBankole     KSPSetType(ksp, KSPPREONLY);
96fcf85c8cSAdelekeBankole     KSPGetPC(ksp, &pc);
97fcf85c8cSAdelekeBankole     PCSetType(pc, PCLU);
98fcf85c8cSAdelekeBankole   }
99fcf85c8cSAdelekeBankole   /*
100fcf85c8cSAdelekeBankole      Set SNES/KSP/KSP/PC runtime options, e.g.,
101fcf85c8cSAdelekeBankole          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
102fcf85c8cSAdelekeBankole      These options will override those specified above as long as
103fcf85c8cSAdelekeBankole      SNESSetFromOptions() is called _after_ any other customization
104fcf85c8cSAdelekeBankole      routines.
105fcf85c8cSAdelekeBankole   */
106fcf85c8cSAdelekeBankole   PetscCall(SNESSetFromOptions(snes));
107fcf85c8cSAdelekeBankole 
108fcf85c8cSAdelekeBankole   PetscCall(SNESSolve(snes, NULL, x));
109fcf85c8cSAdelekeBankole   //PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
110fcf85c8cSAdelekeBankole 
111fcf85c8cSAdelekeBankole   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112fcf85c8cSAdelekeBankole      Free work space.  All PETSc objects should be destroyed when they
113fcf85c8cSAdelekeBankole      are no longer needed.
114fcf85c8cSAdelekeBankole    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115fcf85c8cSAdelekeBankole 
116fcf85c8cSAdelekeBankole   PetscCall(PetscFree2(blasius->x, weight));
117fcf85c8cSAdelekeBankole   PetscCall(PetscFree(blasius));
118fcf85c8cSAdelekeBankole   PetscCall(VecDestroy(&x));
119fcf85c8cSAdelekeBankole   PetscCall(VecDestroy(&r));
120fcf85c8cSAdelekeBankole   PetscCall(SNESDestroy(&snes));
121fcf85c8cSAdelekeBankole   PetscCall(PetscFinalize());
122fcf85c8cSAdelekeBankole   return 0;
123fcf85c8cSAdelekeBankole }
124fcf85c8cSAdelekeBankole 
125*f6dfbefdSBarry Smith /*
126fcf85c8cSAdelekeBankole    Helper function to evaluate Chebyshev polynomials with a set of coefficients
127fcf85c8cSAdelekeBankole    with all their derivatives represented as a recurrence table
128*f6dfbefdSBarry Smith */
129fcf85c8cSAdelekeBankole static void ChebyshevEval(PetscInt N, const PetscScalar *Tf, PetscReal x, PetscReal dx_deta, PetscScalar *f) {
130fcf85c8cSAdelekeBankole   PetscScalar table[4][3] = {
1319371c9d4SSatish Balay     {1, x, 2 * x * x - 1},
1329371c9d4SSatish Balay     {0, 1, 4 * x        },
1339371c9d4SSatish Balay     {0, 0, 4            },
1349371c9d4SSatish Balay     {0, 0, 0            }  /* Chebyshev polynomials T_0, T_1, T_2 of the first kind in (-1,1)  */
135fcf85c8cSAdelekeBankole   };
1369371c9d4SSatish Balay   for (int i = 0; i < 4; i++) { f[i] = table[i][0] * Tf[0] + table[i][1] * Tf[1] + table[i][2] * Tf[2]; /* i-th derivative of f */ }
137fcf85c8cSAdelekeBankole   for (int i = 3; i < N; i++) {
138fcf85c8cSAdelekeBankole     table[0][i % 3] = 2 * x * table[0][(i - 1) % 3] - table[0][(i - 2) % 3]; /* T_n(x) = 2xT_{n-1}(x) - T_{n-2}(x) */
139fcf85c8cSAdelekeBankole     /* Differentiate Chebyshev polynomials with the recurrence relation */
1409371c9d4SSatish Balay     for (int j = 1; j < 4; j++) { table[j][i % 3] = i * (2 * table[j - 1][(i - 1) % 3] + table[j][(i - 2) % 3] / (i - 2)); /* T'_{n}(x)/n = 2T_{n-1}(x) + T'_{n-2}(x)/n-2 */ }
141ad540459SPierre Jolivet     for (int j = 0; j < 4; j++) f[j] += table[j][i % 3] * Tf[i];
142fcf85c8cSAdelekeBankole   }
143fcf85c8cSAdelekeBankole   for (int i = 1; i < 4; i++) {
144fcf85c8cSAdelekeBankole     for (int j = 0; j < i; j++) f[i] *= dx_deta; /* Here happens the physics of the problem */
145fcf85c8cSAdelekeBankole   }
146fcf85c8cSAdelekeBankole }
147fcf85c8cSAdelekeBankole 
148fcf85c8cSAdelekeBankole /*
149fcf85c8cSAdelekeBankole    FormFunction - Evaluates nonlinear function, F(x).
150fcf85c8cSAdelekeBankole 
151fcf85c8cSAdelekeBankole    Input Parameters:
152fcf85c8cSAdelekeBankole .  snes - the SNES context
153fcf85c8cSAdelekeBankole .  X    - input vector
154fcf85c8cSAdelekeBankole .  ctx  - optional user-defined context
155fcf85c8cSAdelekeBankole 
156fcf85c8cSAdelekeBankole    Output Parameter:
157fcf85c8cSAdelekeBankole .  R - function vector
158fcf85c8cSAdelekeBankole  */
1599371c9d4SSatish Balay PetscErrorCode FormFunction(SNES snes, Vec X, Vec R, void *ctx) {
160fcf85c8cSAdelekeBankole   Blasius           *blasius = (Blasius *)ctx;
161fcf85c8cSAdelekeBankole   const PetscScalar *Tf, *Th; /* Tf and Th are Chebyshev coefficients */
162fcf85c8cSAdelekeBankole   PetscScalar       *r, f[4], h[4];
163fcf85c8cSAdelekeBankole   PetscInt           N  = blasius->N;
164fcf85c8cSAdelekeBankole   PetscReal          Ma = blasius->Ma, Pr = blasius->Pr;
165fcf85c8cSAdelekeBankole 
166fcf85c8cSAdelekeBankole   /*
167fcf85c8cSAdelekeBankole    Get pointers to vector data.
168fcf85c8cSAdelekeBankole       - For default PETSc vectors, VecGetArray() returns a pointer to
169fcf85c8cSAdelekeBankole         the data array.  Otherwise, the routine is implementation dependent.
170fcf85c8cSAdelekeBankole       - You MUST call VecRestoreArray() when you no longer need access to
171fcf85c8cSAdelekeBankole         the array.
172fcf85c8cSAdelekeBankole    */
173fcf85c8cSAdelekeBankole   PetscCall(VecGetArrayRead(X, &Tf));
174fcf85c8cSAdelekeBankole   Th = Tf + N;
175fcf85c8cSAdelekeBankole   PetscCall(VecGetArray(R, &r));
176fcf85c8cSAdelekeBankole 
177fcf85c8cSAdelekeBankole   /* Compute function */
178fcf85c8cSAdelekeBankole   ChebyshevEval(N, Tf, -1., blasius->dx_deta, f);
179fcf85c8cSAdelekeBankole   r[0] = f[0];
180fcf85c8cSAdelekeBankole   r[1] = f[1];
181fcf85c8cSAdelekeBankole   ChebyshevEval(N, Tf, 1., blasius->dx_deta, f);
182fcf85c8cSAdelekeBankole   r[2] = f[1] - 1; /* Right end boundary condition */
183fcf85c8cSAdelekeBankole   for (int i = 0; i < N - 3; i++) {
184fcf85c8cSAdelekeBankole     ChebyshevEval(N, Tf, blasius->x[i], blasius->dx_deta, f);
185fcf85c8cSAdelekeBankole     r[3 + i] = 2 * f[3] + f[2] * f[0];
186fcf85c8cSAdelekeBankole     ChebyshevEval(N - 1, Th, blasius->x[i], blasius->dx_deta, h);
187fcf85c8cSAdelekeBankole     r[N + 2 + i] = h[2] + Pr * f[0] * h[1] + Pr * (blasius->gamma - 1) * PetscSqr(Ma * f[2]);
188fcf85c8cSAdelekeBankole   }
189fcf85c8cSAdelekeBankole   ChebyshevEval(N - 1, Th, -1., blasius->dx_deta, h);
190fcf85c8cSAdelekeBankole   r[N] = h[0] - blasius->h_0; /* Left end boundary condition */
191fcf85c8cSAdelekeBankole   ChebyshevEval(N - 1, Th, 1., blasius->dx_deta, h);
192fcf85c8cSAdelekeBankole   r[N + 1] = h[0] - 1; /* Left end boundary condition */
193fcf85c8cSAdelekeBankole 
194fcf85c8cSAdelekeBankole   /* Restore vectors */
195fcf85c8cSAdelekeBankole   PetscCall(VecRestoreArrayRead(X, &Tf));
196fcf85c8cSAdelekeBankole   PetscCall(VecRestoreArray(R, &r));
197fcf85c8cSAdelekeBankole   return 0;
198fcf85c8cSAdelekeBankole }
199fcf85c8cSAdelekeBankole 
200fcf85c8cSAdelekeBankole /*TEST
201fcf85c8cSAdelekeBankole 
202fcf85c8cSAdelekeBankole    test:
203fcf85c8cSAdelekeBankole       args: -snes_monitor -pc_type svd
204fcf85c8cSAdelekeBankole       requires: !single
205fcf85c8cSAdelekeBankole 
206fcf85c8cSAdelekeBankole TEST*/
207