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