11ef27442SStefano Zampini #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 2b79b07cfSJed Brown 31ef27442SStefano Zampini typedef struct { 41ef27442SStefano Zampini PetscBool transpose_solve; 51ef27442SStefano Zampini } SNES_KSPONLY; 6b79b07cfSJed Brown 7d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_KSPONLY(SNES snes) 8d71ae5a4SJacob Faibussowitsch { 91ef27442SStefano Zampini SNES_KSPONLY *ksponly = (SNES_KSPONLY *)snes->data; 10b79b07cfSJed Brown PetscInt lits; 11b79b07cfSJed Brown Vec Y, X, F; 12b79b07cfSJed Brown 13b79b07cfSJed Brown PetscFunctionBegin; 140b121fc5SBarry Smith PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 15c579b300SPatrick Farrell 16b79b07cfSJed Brown snes->numFailures = 0; 17b79b07cfSJed Brown snes->numLinearSolveFailures = 0; 18b79b07cfSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 19b79b07cfSJed Brown snes->iter = 0; 20b79b07cfSJed Brown snes->norm = 0.0; 21b79b07cfSJed Brown 22b79b07cfSJed Brown X = snes->vec_sol; 23b79b07cfSJed Brown F = snes->vec_func; 24b79b07cfSJed Brown Y = snes->vec_sol_update; 25b79b07cfSJed Brown 26451891c6SStefano Zampini if (!snes->vec_func_init_set) { 279566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 28451891c6SStefano Zampini } else snes->vec_func_init_set = PETSC_FALSE; 29451891c6SStefano Zampini 30a5eaddd9SBarry Smith if (snes->numbermonitors) { 31a5eaddd9SBarry Smith PetscReal fnorm; 329566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 339566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm)); 34a5eaddd9SBarry Smith } 35b79b07cfSJed Brown 365341784dSBarry Smith /* Call general purpose update function */ 37dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, 0); 385341784dSBarry Smith 39b79b07cfSJed Brown /* Solve J Y = F, where J is Jacobian matrix */ 409566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 4107b62357SFande Kong 4207b62357SFande Kong SNESCheckJacobianDomainerror(snes); 4307b62357SFande Kong 449566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 451ef27442SStefano Zampini if (ksponly->transpose_solve) { 469566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(snes->ksp, F, Y)); 471ef27442SStefano Zampini } else { 489566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, Y)); 491ef27442SStefano Zampini } 50422a814eSBarry Smith snes->reason = SNES_CONVERGED_ITS; 51422a814eSBarry Smith SNESCheckKSPSolve(snes); 521aa26658SKarl Rupp 539566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 5463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 55b79b07cfSJed Brown snes->iter++; 56b79b07cfSJed Brown 57b79b07cfSJed Brown /* Take the computed step. */ 589566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, -1.0, Y)); 59a5eaddd9SBarry Smith if (snes->numbermonitors) { 60a5eaddd9SBarry Smith PetscReal fnorm; 619566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 629566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 639566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 1, fnorm)); 64a5eaddd9SBarry Smith } 65*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 66b79b07cfSJed Brown } 67b79b07cfSJed Brown 68d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_KSPONLY(SNES snes) 69d71ae5a4SJacob Faibussowitsch { 70b79b07cfSJed Brown PetscFunctionBegin; 719566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 72*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 73b79b07cfSJed Brown } 74b79b07cfSJed Brown 75d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_KSPONLY(SNES snes) 76d71ae5a4SJacob Faibussowitsch { 77b79b07cfSJed Brown PetscFunctionBegin; 789566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 79*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 80b79b07cfSJed Brown } 81b79b07cfSJed Brown 82b79b07cfSJed Brown /*MC 83f6dfbefdSBarry Smith SNESKSPONLY - Nonlinear solver that performs one Newton step and does not compute any norms. 84f6dfbefdSBarry Smith The main purpose of this solver is to solve linear problems using the `SNES` interface, without 85b79b07cfSJed Brown any additional overhead in the form of vector operations. 86b79b07cfSJed Brown 87b79b07cfSJed Brown Level: beginner 88b79b07cfSJed Brown 89f6dfbefdSBarry Smith .seealso: `SNES`, `SNESType`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR` 90b79b07cfSJed Brown M*/ 91d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes) 92d71ae5a4SJacob Faibussowitsch { 931ef27442SStefano Zampini SNES_KSPONLY *ksponly; 94b79b07cfSJed Brown 95b79b07cfSJed Brown PetscFunctionBegin; 96b79b07cfSJed Brown snes->ops->setup = SNESSetUp_KSPONLY; 97b79b07cfSJed Brown snes->ops->solve = SNESSolve_KSPONLY; 98b79b07cfSJed Brown snes->ops->destroy = SNESDestroy_KSPONLY; 999e5d0892SLisandro Dalcin snes->ops->setfromoptions = NULL; 1009e5d0892SLisandro Dalcin snes->ops->view = NULL; 1019e5d0892SLisandro Dalcin snes->ops->reset = NULL; 102b79b07cfSJed Brown 10342f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 104efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 10542f4f86dSBarry Smith 1064fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 1074fc747eaSLawrence Mitchell 1084dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ksponly)); 1091ef27442SStefano Zampini snes->data = (void *)ksponly; 110*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1111ef27442SStefano Zampini } 1121ef27442SStefano Zampini 1131ef27442SStefano Zampini /*MC 114f6dfbefdSBarry Smith SNESKSPTRANSPOSEONLY - Nonlinear solver that performs one Newton step and does not compute any norms. 115f6dfbefdSBarry Smith The main purpose of this solver is to solve transposed linear problems using the `SNES` interface, without 1161ef27442SStefano Zampini any additional overhead in the form of vector operations within adjoint solvers. 1171ef27442SStefano Zampini 1181ef27442SStefano Zampini Level: beginner 1191ef27442SStefano Zampini 120f6dfbefdSBarry Smith .seealso: `SNES`, `SNESType`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESKSPTRANSPOSEONLY`, `SNESNEWTONLS`, `SNESNEWTONTR` 1211ef27442SStefano Zampini M*/ 122d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_KSPTRANSPOSEONLY(SNES snes) 123d71ae5a4SJacob Faibussowitsch { 1241ef27442SStefano Zampini SNES_KSPONLY *kspo; 1251ef27442SStefano Zampini 1261ef27442SStefano Zampini PetscFunctionBegin; 1279566063dSJacob Faibussowitsch PetscCall(SNESCreate_KSPONLY(snes)); 1289566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)snes, SNESKSPTRANSPOSEONLY)); 1291ef27442SStefano Zampini kspo = (SNES_KSPONLY *)snes->data; 1301ef27442SStefano Zampini kspo->transpose_solve = PETSC_TRUE; 131*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 132b79b07cfSJed Brown } 133