source: CIVL/examples/compare/PETSc/ex2d.c

main
Last change on this file 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: 3.0 KB
RevLine 
[1f51041]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
[28be836]17#include "petsc.h"
[1f51041]18
19typedef struct {
20 PassiveReal param; /* test problem parameter */
21} AppCtx;
22
[28be836]23PetscScalar FormFunctionNonlin(PetscScalar ,PetscScalar ,PetscScalar ,PetscScalar ,PetscScalar ,PetscReal ,PetscReal ,PetscReal );
24PetscScalar FormFunctionLin(PetscScalar ,PetscScalar ,PetscScalar ,PetscScalar ,PetscScalar ,PetscReal ,PetscReal ,PetscReal );
25PetscScalar FormFunctionPt(PetscScalar ,PetscScalar ,PetscScalar ,PetscScalar ,PetscScalar ,PetscReal ,PetscReal ,PetscReal );
[1f51041]26/* ------------------------------------------------------------------- */
27#undef __FUNCT__
28#define __FUNCT__ "FormFunctionLocal"
29/*
30 FormFunctionLocal - Evaluates nonlinear function, F(x) on local process patch
31 */
32PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
33{
34 PetscErrorCode ierr;
35 PetscInt i,j;
36 PetscReal lambda,hx,hy,hxdhy,hydhx,sc;
37 PetscScalar u,uxx,uyy;
38
39 PetscFunctionBegin;
40
41 lambda = user->param;
42 hx = 1.0/(PetscReal)(info->mx-1);
43 hy = 1.0/(PetscReal)(info->my-1);
44 sc = hx*hy*lambda;
45 hxdhy = hx/hy;
46 hydhx = hy/hx;
47
48 /*
49 Compute function over the locally owned part of the grid
50 */
51 for (j=info->ys; j<info->ys+info->ym; j++) {
52 for (i=info->xs; i<info->xs+info->xm; i++) {
53 if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
54 f[j][i] = 2.0*(hydhx+hxdhy)*x[j][i];
55 ierr = PetscLogFlops(3.0);CHKERRQ(ierr);
56 } else {
[28be836]57 f[j][i] = FormFunctionPt(x[j][i],x[j-1][i],x[j+1][i],x[j][i-1],x[j][i+1],hxdhy,hydhx,sc);
[1f51041]58 ierr = PetscLogFlops(11.0);CHKERRQ(ierr);
59 }
60 }
61 }
62 PetscFunctionReturn(0);
63}
64
65PetscScalar FormFunctionPt(PetscScalar C,PetscScalar S,PetscScalar N,PetscScalar W,PetscScalar E,PetscReal hxdhy,PetscReal hydhx,PetscReal sc) {
66 PetscScalar f;
67
[28be836]68 f = FormFunctionLin(C,S,N,W,E,hxdhy,hydhx,sc);
[3073847]69 f -= FormFunctionNonlin(C,S,N,W,E,hxdhy,hydhx,sc);
[1f51041]70
71 return f;
72}
73
74PetscScalar FormFunctionLin(PetscScalar C,PetscScalar S,PetscScalar N,PetscScalar W,PetscScalar E,PetscReal hxdhy,PetscReal hydhx,PetscReal sc) {
75 PetscScalar uxx,uyy;
76
77 uxx = (- W - E)*hydhx;
78 uyy = (- S - N)*hxdhy;
79 return (uxx + uyy);
80}
81
82PetscScalar FormFunctionNonlin(PetscScalar C,PetscScalar S,PetscScalar N,PetscScalar W,PetscScalar E,PetscReal hxdhy,PetscReal hydhx,PetscReal sc) {
83 PetscScalar uxx,uyy;
84
85 uxx = (2.0*C)*hydhx;
86 uyy = (2.0*C)*hxdhy;
87 return (uxx + uyy + sc*PetscExpScalar(C));
88}
Note: See TracBrowser for help on using the repository browser.