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 PetscInt its; 25c4762a1bSJed Brown AppCtx user; 26c4762a1bSJed Brown PetscMPIInt size; 27c4762a1bSJed Brown 289566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 30*08401ef6SPierre Jolivet PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* Allocate vectors / matrix */ 339566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&x)); 349566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x,PETSC_DECIDE,1)); 359566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&r)); 37c4762a1bSJed Brown 389566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD,1,1,1,NULL,&J)); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Create / set-up SNES */ 419566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes)); 429566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,r,UserFunction,&user)); 439566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,J,J,UserJacobian,&user)); 449566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* Set initial guess (=1) and target value */ 47c4762a1bSJed Brown user.value = 1e-4; 48c4762a1bSJed Brown 499566063dSJacob Faibussowitsch PetscCall(VecSet(x,1.0)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Set initial guess / solve */ 529566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes,NULL,x)); 539566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&its)); 549566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its)); 559566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* Done */ 589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 619566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 629566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 63b122ec5aSJacob Faibussowitsch return 0; 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* 67c4762a1bSJed Brown UserFunction - for nonlinear function x^2 - value = 0 68c4762a1bSJed Brown */ 69c4762a1bSJed Brown PetscErrorCode UserFunction(SNES snes,Vec X,Vec F,void *ptr) 70c4762a1bSJed Brown { 71c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 72c4762a1bSJed Brown PetscInt N,i; 73c4762a1bSJed Brown PetscScalar *f; 74c4762a1bSJed Brown PetscReal half; 75c4762a1bSJed Brown const PetscScalar *x; 76c4762a1bSJed Brown 77c4762a1bSJed Brown half = 0.5; 78c4762a1bSJed Brown 799566063dSJacob Faibussowitsch PetscCall(VecGetSize(X,&N)); 809566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 819566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* Calculate residual */ 84c4762a1bSJed Brown for (i=0; i<N; ++i) { 85c4762a1bSJed Brown /* 86c4762a1bSJed Brown Test for domain error. 87c4762a1bSJed Brown Artifical test is applied. With starting value 1.0, first iterate will be 0.5 + user->value/2. 88c4762a1bSJed Brown Declare (0.5-value,0.5+value) to be infeasible. 89c4762a1bSJed Brown In later iterations, snes->domainerror should be cleared, allowing iterations in the feasible region to be accepted. 90c4762a1bSJed Brown */ 91c4762a1bSJed Brown if ((half-user->value) < PetscRealPart(x[i]) && PetscRealPart(x[i]) < (half+user->value)) { 929566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"DOMAIN ERROR: x=%g\n",(double)PetscRealPart(x[i]))); 939566063dSJacob Faibussowitsch PetscCall(SNESSetFunctionDomainError(snes)); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown f[i] = x[i]*x[i] - user->value; 96c4762a1bSJed Brown } 979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 989566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 99c4762a1bSJed Brown return 0; 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* 103c4762a1bSJed Brown UserJacobian - for nonlinear function x^2 - value = 0 104c4762a1bSJed Brown */ 105c4762a1bSJed Brown PetscErrorCode UserJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr) 106c4762a1bSJed Brown { 107c4762a1bSJed Brown PetscInt N,i,row,col; 108c4762a1bSJed Brown const PetscScalar *x; 109c4762a1bSJed Brown PetscScalar v; 110c4762a1bSJed Brown 1119566063dSJacob Faibussowitsch PetscCall(VecGetSize(X,&N)); 1129566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* Calculate Jacobian */ 115c4762a1bSJed Brown for (i=0; i<N; ++i) { 116c4762a1bSJed Brown row = i; 117c4762a1bSJed Brown col = i; 118c4762a1bSJed Brown v = 2*x[i]; 1199566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac,1,&row,1,&col,&v,INSERT_VALUES)); 120c4762a1bSJed Brown } 1219566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 1229566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 1239566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown if (jac != J) { 1269566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 1279566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 128c4762a1bSJed Brown } 129c4762a1bSJed Brown return 0; 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132c4762a1bSJed Brown /*TEST 133c4762a1bSJed Brown 134c4762a1bSJed Brown build: 135dfd57a17SPierre Jolivet requires: !single !defined(PETSC_HAVE_SUN_CXX) !complex 136c4762a1bSJed Brown 137c4762a1bSJed Brown test: 138c4762a1bSJed Brown args: -snes_monitor_solution -snes_linesearch_monitor 139c4762a1bSJed Brown 140c4762a1bSJed Brown TEST*/ 141