source: CIVL/examples/compare/PETSc/ex2a.c@ bb03188

main test-branch
Last change on this file since bb03188 was ea777aa, checked in by Alex Wilton <awilton@…>, 3 years ago

Moved examples, include, build_default.properties, common.xml, and README out from dev.civl.com into the root of the repo.

git-svn-id: svn://vsl.cis.udel.edu/civl/trunk@5704 fb995dde-84ed-4084-dfe6-e5aef3e2452c

  • Property mode set to 100644
File size: 1.9 KB
Line 
1/* ------------------------------------------------------------------------
2
3 Solid Fuel Ignition (SFI) problem. This problem is modeled by
4 the partial differential equation
5
6 -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
7
8 with boundary conditions
9
10 u = 0 for x = 0, x = 1, y = 0, y = 1.
11
12 A finite difference approximation with the usual 5-point stencil
13 is used to discretize the boundary value problem to obtain a nonlinear
14 system of equations.
15 ------------------------------------------------------------------------- */
16
17#include "petsc.h"
18
19typedef struct {
20 PassiveReal param; /* test problem parameter */
21} AppCtx;
22
23/* ------------------------------------------------------------------- */
24#undef __FUNCT__
25#define __FUNCT__ "FormFunctionLocal"
26/*
27 FormFunctionLocal - Evaluates nonlinear function, F(x) on local process patch
28 */
29PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
30{
31 PetscErrorCode ierr;
32 PetscInt i,j;
33 PetscReal lambda,hx,hy,hxdhy,hydhx,sc;
34 PetscScalar u,uxx,uyy;
35
36 PetscFunctionBegin;
37
38 lambda = user->param;
39 hx = 1.0/(PetscReal)(info->mx-1);
40 hy = 1.0/(PetscReal)(info->my-1);
41 sc = hx*hy*lambda;
42 hxdhy = hx/hy;
43 hydhx = hy/hx;
44 /*
45 Compute function over the locally owned part of the grid
46 */
47 for (j=info->ys; j<info->ys+info->ym; j++) {
48 for (i=info->xs; i<info->xs+info->xm; i++) {
49 if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
50 f[j][i] = 2.0*(hydhx+hxdhy)*x[j][i];
51 } else {
52 u = x[j][i];
53 uxx = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx;
54 uyy = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy;
55 f[j][i] = uxx + uyy - sc*PetscExpScalar(u);
56 }
57 }
58 }
59 ierr = PetscLogFlops(11.0*info->ym*info->xm);CHKERRQ(ierr);
60 PetscFunctionReturn(0);
61}
Note: See TracBrowser for help on using the repository browser.