1*c4762a1bSJed Brown static char help[] = "Artificial test to check that snes->domainerror is being reset appropriately"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown /* ------------------------------------------------------------------------ 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown Artificial test to check that snes->domainerror is being reset appropriately 6*c4762a1bSJed Brown 7*c4762a1bSJed Brown ------------------------------------------------------------------------- */ 8*c4762a1bSJed Brown 9*c4762a1bSJed Brown #define PETSC_SKIP_COMPLEX 10*c4762a1bSJed Brown #include <petscsnes.h> 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown typedef struct { 13*c4762a1bSJed Brown PetscReal value; /* parameter in nonlinear function */ 14*c4762a1bSJed Brown } AppCtx; 15*c4762a1bSJed Brown 16*c4762a1bSJed Brown PetscErrorCode UserFunction(SNES,Vec,Vec,void*); 17*c4762a1bSJed Brown PetscErrorCode UserJacobian(SNES,Vec,Mat,Mat,void*); 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown int main(int argc,char **argv) 20*c4762a1bSJed Brown { 21*c4762a1bSJed Brown SNES snes; 22*c4762a1bSJed Brown Vec x,r; 23*c4762a1bSJed Brown Mat J; 24*c4762a1bSJed Brown PetscErrorCode ierr; 25*c4762a1bSJed Brown PetscInt its; 26*c4762a1bSJed Brown AppCtx user; 27*c4762a1bSJed Brown PetscMPIInt size; 28*c4762a1bSJed Brown 29*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 30*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 31*c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown /* Allocate vectors / matrix */ 34*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 35*c4762a1bSJed Brown ierr = VecSetSizes(x,PETSC_DECIDE,1);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = VecSetFromOptions(x);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 38*c4762a1bSJed Brown 39*c4762a1bSJed Brown ierr = MatCreateSeqAIJ(PETSC_COMM_WORLD,1,1,1,NULL,&J);CHKERRQ(ierr); 40*c4762a1bSJed Brown 41*c4762a1bSJed Brown /* Create / set-up SNES */ 42*c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 43*c4762a1bSJed Brown ierr = SNESSetFunction(snes,r,UserFunction,&user);CHKERRQ(ierr); 44*c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,UserJacobian,&user);CHKERRQ(ierr); 45*c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown /* Set initial guess (=1) and target value */ 48*c4762a1bSJed Brown user.value = 1e-4; 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown ierr = VecSet(x,1.0);CHKERRQ(ierr); 51*c4762a1bSJed Brown 52*c4762a1bSJed Brown /* Set initial guess / solve */ 53*c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 57*c4762a1bSJed Brown 58*c4762a1bSJed Brown /* Done */ 59*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = PetscFinalize(); 64*c4762a1bSJed Brown return ierr; 65*c4762a1bSJed Brown } 66*c4762a1bSJed Brown 67*c4762a1bSJed Brown /* 68*c4762a1bSJed Brown UserFunction - for nonlinear function x^2 - value = 0 69*c4762a1bSJed Brown */ 70*c4762a1bSJed Brown PetscErrorCode UserFunction(SNES snes,Vec X,Vec F,void *ptr) 71*c4762a1bSJed Brown { 72*c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 73*c4762a1bSJed Brown PetscInt N,i; 74*c4762a1bSJed Brown PetscScalar *f; 75*c4762a1bSJed Brown PetscReal half; 76*c4762a1bSJed Brown const PetscScalar *x; 77*c4762a1bSJed Brown PetscErrorCode ierr; 78*c4762a1bSJed Brown 79*c4762a1bSJed Brown half = 0.5; 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown ierr = VecGetSize(X,&N);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 84*c4762a1bSJed Brown 85*c4762a1bSJed Brown /* Calculate residual */ 86*c4762a1bSJed Brown for(i=0; i<N; ++i) { 87*c4762a1bSJed Brown /* 88*c4762a1bSJed Brown Test for domain error. 89*c4762a1bSJed Brown Artifical test is applied. With starting value 1.0, first iterate will be 0.5 + user->value/2. 90*c4762a1bSJed Brown Declare (0.5-value,0.5+value) to be infeasible. 91*c4762a1bSJed Brown In later iterations, snes->domainerror should be cleared, allowing iterations in the feasible region to be accepted. 92*c4762a1bSJed Brown */ 93*c4762a1bSJed Brown if( (half-user->value) < PetscRealPart(x[i]) && PetscRealPart(x[i]) < (half+user->value) ) { 94*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"DOMAIN ERROR: x=%g\n",(double)PetscRealPart(x[i]));CHKERRQ(ierr); 95*c4762a1bSJed Brown ierr = SNESSetFunctionDomainError(snes);CHKERRQ(ierr); 96*c4762a1bSJed Brown } 97*c4762a1bSJed Brown f[i] = x[i]*x[i] - user->value; 98*c4762a1bSJed Brown } 99*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 100*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 101*c4762a1bSJed Brown return 0; 102*c4762a1bSJed Brown } 103*c4762a1bSJed Brown 104*c4762a1bSJed Brown /* 105*c4762a1bSJed Brown UserJacobian - for nonlinear function x^2 - value = 0 106*c4762a1bSJed Brown */ 107*c4762a1bSJed Brown PetscErrorCode UserJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr) 108*c4762a1bSJed Brown { 109*c4762a1bSJed Brown PetscInt N,i,row,col; 110*c4762a1bSJed Brown const PetscScalar *x; 111*c4762a1bSJed Brown PetscScalar v; 112*c4762a1bSJed Brown PetscErrorCode ierr; 113*c4762a1bSJed Brown 114*c4762a1bSJed Brown ierr = VecGetSize(X,&N);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown /* Calculate Jacobian */ 118*c4762a1bSJed Brown for (i=0; i<N; ++i) { 119*c4762a1bSJed Brown row = i; 120*c4762a1bSJed Brown col = i; 121*c4762a1bSJed Brown v = 2*x[i]; 122*c4762a1bSJed Brown ierr = MatSetValues(jac,1,&row,1,&col,&v,INSERT_VALUES);CHKERRQ(ierr); 123*c4762a1bSJed Brown } 124*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 125*c4762a1bSJed Brown ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 127*c4762a1bSJed Brown 128*c4762a1bSJed Brown if (jac != J) { 129*c4762a1bSJed Brown ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 130*c4762a1bSJed Brown ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 131*c4762a1bSJed Brown } 132*c4762a1bSJed Brown return 0; 133*c4762a1bSJed Brown } 134*c4762a1bSJed Brown 135*c4762a1bSJed Brown /*TEST 136*c4762a1bSJed Brown 137*c4762a1bSJed Brown build: 138*c4762a1bSJed Brown requires: !single !define(PETSC_HAVE_SUN_CXX) !complex 139*c4762a1bSJed Brown 140*c4762a1bSJed Brown test: 141*c4762a1bSJed Brown args: -snes_monitor_solution -snes_linesearch_monitor 142*c4762a1bSJed Brown 143*c4762a1bSJed Brown TEST*/ 144*c4762a1bSJed Brown 145