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 9*74c01ffaSSatish Balay #if !defined(PETSC_SKIP_COMPLEX) 10c4762a1bSJed Brown #define PETSC_SKIP_COMPLEX 11*74c01ffaSSatish Balay #endif 12c4762a1bSJed Brown #include <petscsnes.h> 13c4762a1bSJed Brown 14c4762a1bSJed Brown typedef struct { 15c4762a1bSJed Brown PetscReal value; /* parameter in nonlinear function */ 16c4762a1bSJed Brown } AppCtx; 17c4762a1bSJed Brown 18c4762a1bSJed Brown PetscErrorCode UserFunction(SNES, Vec, Vec, void *); 19c4762a1bSJed Brown PetscErrorCode UserJacobian(SNES, Vec, Mat, Mat, void *); 20c4762a1bSJed Brown 21d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 22d71ae5a4SJacob Faibussowitsch { 23c4762a1bSJed Brown SNES snes; 24c4762a1bSJed Brown Vec x, r; 25c4762a1bSJed Brown Mat J; 26c4762a1bSJed Brown PetscInt its; 27c4762a1bSJed Brown AppCtx user; 28c4762a1bSJed Brown PetscMPIInt size; 29c4762a1bSJed Brown 30327415f7SBarry Smith PetscFunctionBeginUser; 319566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 3308401ef6SPierre Jolivet PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 34c4762a1bSJed Brown 35c4762a1bSJed Brown /* Allocate vectors / matrix */ 369566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 379566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, PETSC_DECIDE, 1)); 389566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 399566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 40c4762a1bSJed Brown 419566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, 1, 1, 1, NULL, &J)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Create / set-up SNES */ 449566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 459566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, UserFunction, &user)); 469566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, UserJacobian, &user)); 479566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* Set initial guess (=1) and target value */ 50c4762a1bSJed Brown user.value = 1e-4; 51c4762a1bSJed Brown 529566063dSJacob Faibussowitsch PetscCall(VecSet(x, 1.0)); 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* Set initial guess / solve */ 559566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x)); 569566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 5763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its)); 589566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Done */ 619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 639566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 649566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 659566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 66b122ec5aSJacob Faibussowitsch return 0; 67c4762a1bSJed Brown } 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* 70c4762a1bSJed Brown UserFunction - for nonlinear function x^2 - value = 0 71c4762a1bSJed Brown */ 72d71ae5a4SJacob Faibussowitsch PetscErrorCode UserFunction(SNES snes, Vec X, Vec F, void *ptr) 73d71ae5a4SJacob Faibussowitsch { 74c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 75c4762a1bSJed Brown PetscInt N, i; 76c4762a1bSJed Brown PetscScalar *f; 77c4762a1bSJed Brown PetscReal half; 78c4762a1bSJed Brown const PetscScalar *x; 79c4762a1bSJed Brown 80c4762a1bSJed Brown half = 0.5; 813ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 829566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &N)); 839566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 849566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Calculate residual */ 87c4762a1bSJed Brown for (i = 0; i < N; ++i) { 88c4762a1bSJed Brown /* 89c4762a1bSJed Brown Test for domain error. 90d5b43468SJose E. Roman Artificial test is applied. With starting value 1.0, first iterate will be 0.5 + user->value/2. 91c4762a1bSJed Brown Declare (0.5-value,0.5+value) to be infeasible. 92c4762a1bSJed Brown In later iterations, snes->domainerror should be cleared, allowing iterations in the feasible region to be accepted. 93c4762a1bSJed Brown */ 94c4762a1bSJed Brown if ((half - user->value) < PetscRealPart(x[i]) && PetscRealPart(x[i]) < (half + user->value)) { 959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "DOMAIN ERROR: x=%g\n", (double)PetscRealPart(x[i]))); 969566063dSJacob Faibussowitsch PetscCall(SNESSetFunctionDomainError(snes)); 97c4762a1bSJed Brown } 98c4762a1bSJed Brown f[i] = x[i] * x[i] - user->value; 99c4762a1bSJed Brown } 1009566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 1019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 1023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* 106c4762a1bSJed Brown UserJacobian - for nonlinear function x^2 - value = 0 107c4762a1bSJed Brown */ 108d71ae5a4SJacob Faibussowitsch PetscErrorCode UserJacobian(SNES snes, Vec X, Mat J, Mat jac, void *ptr) 109d71ae5a4SJacob Faibussowitsch { 110c4762a1bSJed Brown PetscInt N, i, row, col; 111c4762a1bSJed Brown const PetscScalar *x; 112c4762a1bSJed Brown PetscScalar v; 113c4762a1bSJed Brown 1143ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1159566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &N)); 1169566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* Calculate Jacobian */ 119c4762a1bSJed Brown for (i = 0; i < N; ++i) { 120c4762a1bSJed Brown row = i; 121c4762a1bSJed Brown col = i; 122c4762a1bSJed Brown v = 2 * x[i]; 1239566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 1, &row, 1, &col, &v, INSERT_VALUES)); 124c4762a1bSJed Brown } 1259566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 1269566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 1279566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 128c4762a1bSJed Brown 129c4762a1bSJed Brown if (jac != J) { 1309566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 1319566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 132c4762a1bSJed Brown } 1333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 134c4762a1bSJed Brown } 135c4762a1bSJed Brown 136c4762a1bSJed Brown /*TEST 137c4762a1bSJed Brown 138c4762a1bSJed Brown build: 139dfd57a17SPierre Jolivet requires: !single !defined(PETSC_HAVE_SUN_CXX) !complex 140c4762a1bSJed Brown 141c4762a1bSJed Brown test: 142c4762a1bSJed Brown args: -snes_monitor_solution -snes_linesearch_monitor 143c4762a1bSJed Brown 144c4762a1bSJed Brown TEST*/ 145