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 7b79b07cfSJed Brown static PetscErrorCode SNESSolve_KSPONLY(SNES snes) 8b79b07cfSJed Brown { 91ef27442SStefano Zampini SNES_KSPONLY *ksponly = (SNES_KSPONLY*)snes->data; 10b79b07cfSJed Brown PetscErrorCode ierr; 11b79b07cfSJed Brown PetscInt lits; 12b79b07cfSJed Brown Vec Y,X,F; 13b79b07cfSJed Brown 14b79b07cfSJed Brown PetscFunctionBegin; 15*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(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); 16c579b300SPatrick Farrell 17b79b07cfSJed Brown snes->numFailures = 0; 18b79b07cfSJed Brown snes->numLinearSolveFailures = 0; 19b79b07cfSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 20b79b07cfSJed Brown snes->iter = 0; 21b79b07cfSJed Brown snes->norm = 0.0; 22b79b07cfSJed Brown 23b79b07cfSJed Brown X = snes->vec_sol; 24b79b07cfSJed Brown F = snes->vec_func; 25b79b07cfSJed Brown Y = snes->vec_sol_update; 26b79b07cfSJed Brown 27451891c6SStefano Zampini if (!snes->vec_func_init_set) { 28b79b07cfSJed Brown ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 29451891c6SStefano Zampini } else snes->vec_func_init_set = PETSC_FALSE; 30451891c6SStefano Zampini 31a5eaddd9SBarry Smith if (snes->numbermonitors) { 32a5eaddd9SBarry Smith PetscReal fnorm; 33a5eaddd9SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 34a5eaddd9SBarry Smith ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 35a5eaddd9SBarry Smith } 36b79b07cfSJed Brown 375341784dSBarry Smith /* Call general purpose update function */ 385341784dSBarry Smith if (snes->ops->update) { 395341784dSBarry Smith ierr = (*snes->ops->update)(snes, 0);CHKERRQ(ierr); 405341784dSBarry Smith } 415341784dSBarry Smith 42b79b07cfSJed Brown /* Solve J Y = F, where J is Jacobian matrix */ 43d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 4407b62357SFande Kong 4507b62357SFande Kong SNESCheckJacobianDomainerror(snes); 4607b62357SFande Kong 4723ee1639SBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 481ef27442SStefano Zampini if (ksponly->transpose_solve) { 491ef27442SStefano Zampini ierr = KSPSolveTranspose(snes->ksp,F,Y);CHKERRQ(ierr); 501ef27442SStefano Zampini } else { 51b79b07cfSJed Brown ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr); 521ef27442SStefano Zampini } 53422a814eSBarry Smith snes->reason = SNES_CONVERGED_ITS; 54422a814eSBarry Smith SNESCheckKSPSolve(snes); 551aa26658SKarl Rupp 56b79b07cfSJed Brown ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 577d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 58b79b07cfSJed Brown snes->iter++; 59b79b07cfSJed Brown 60b79b07cfSJed Brown /* Take the computed step. */ 61b79b07cfSJed Brown ierr = VecAXPY(X,-1.0,Y);CHKERRQ(ierr); 62a5eaddd9SBarry Smith if (snes->numbermonitors) { 63a5eaddd9SBarry Smith PetscReal fnorm; 64fef909d3SLawrence Mitchell ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 65a5eaddd9SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 66a5eaddd9SBarry Smith ierr = SNESMonitor(snes,1,fnorm);CHKERRQ(ierr); 67a5eaddd9SBarry Smith } 68b79b07cfSJed Brown PetscFunctionReturn(0); 69b79b07cfSJed Brown } 70b79b07cfSJed Brown 71b79b07cfSJed Brown static PetscErrorCode SNESSetUp_KSPONLY(SNES snes) 72b79b07cfSJed Brown { 736cab3a1bSJed Brown PetscErrorCode ierr; 74b79b07cfSJed Brown 75b79b07cfSJed Brown PetscFunctionBegin; 766cab3a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 77b79b07cfSJed Brown PetscFunctionReturn(0); 78b79b07cfSJed Brown } 79b79b07cfSJed Brown 80b79b07cfSJed Brown static PetscErrorCode SNESDestroy_KSPONLY(SNES snes) 81b79b07cfSJed Brown { 821ef27442SStefano Zampini PetscErrorCode ierr; 83b79b07cfSJed Brown 84b79b07cfSJed Brown PetscFunctionBegin; 851ef27442SStefano Zampini ierr = PetscFree(snes->data);CHKERRQ(ierr); 86b79b07cfSJed Brown PetscFunctionReturn(0); 87b79b07cfSJed Brown } 88b79b07cfSJed Brown 89b79b07cfSJed Brown /* -------------------------------------------------------------------------- */ 90b79b07cfSJed Brown /*MC 91b79b07cfSJed Brown SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms. 92b79b07cfSJed Brown The main purpose of this solver is to solve linear problems using the SNES interface, without 93b79b07cfSJed Brown any additional overhead in the form of vector operations. 94b79b07cfSJed Brown 95b79b07cfSJed Brown Level: beginner 96b79b07cfSJed Brown 9704d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR 98b79b07cfSJed Brown M*/ 998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes) 100b79b07cfSJed Brown { 1011ef27442SStefano Zampini SNES_KSPONLY *ksponly; 1021ef27442SStefano Zampini PetscErrorCode ierr; 103b79b07cfSJed Brown 104b79b07cfSJed Brown PetscFunctionBegin; 105b79b07cfSJed Brown snes->ops->setup = SNESSetUp_KSPONLY; 106b79b07cfSJed Brown snes->ops->solve = SNESSolve_KSPONLY; 107b79b07cfSJed Brown snes->ops->destroy = SNESDestroy_KSPONLY; 1089e5d0892SLisandro Dalcin snes->ops->setfromoptions = NULL; 1099e5d0892SLisandro Dalcin snes->ops->view = NULL; 1109e5d0892SLisandro Dalcin snes->ops->reset = NULL; 111b79b07cfSJed Brown 11242f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 113efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 11442f4f86dSBarry Smith 1154fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 1164fc747eaSLawrence Mitchell 1171ef27442SStefano Zampini ierr = PetscNewLog(snes,&ksponly);CHKERRQ(ierr); 1181ef27442SStefano Zampini snes->data = (void*)ksponly; 1191ef27442SStefano Zampini PetscFunctionReturn(0); 1201ef27442SStefano Zampini } 1211ef27442SStefano Zampini 1221ef27442SStefano Zampini /*MC 1231ef27442SStefano Zampini SNESKSPTRANSPOSEONLY - Nonlinear solver that only performs one Newton step and does not compute any norms. 1241ef27442SStefano Zampini The main purpose of this solver is to solve transposed linear problems using the SNES interface, without 1251ef27442SStefano Zampini any additional overhead in the form of vector operations within adjoint solvers. 1261ef27442SStefano Zampini 1271ef27442SStefano Zampini Level: beginner 1281ef27442SStefano Zampini 1291ef27442SStefano Zampini .seealso: SNESCreate(), SNES, SNESSetType(), SNESKSPTRANSPOSEONLY, SNESNEWTONLS, SNESNEWTONTR 1301ef27442SStefano Zampini M*/ 1311ef27442SStefano Zampini PETSC_EXTERN PetscErrorCode SNESCreate_KSPTRANSPOSEONLY(SNES snes) 1321ef27442SStefano Zampini { 1331ef27442SStefano Zampini SNES_KSPONLY *kspo; 1341ef27442SStefano Zampini PetscErrorCode ierr; 1351ef27442SStefano Zampini 1361ef27442SStefano Zampini PetscFunctionBegin; 1371ef27442SStefano Zampini ierr = SNESCreate_KSPONLY(snes);CHKERRQ(ierr); 13807d8e96eSStefano Zampini ierr = PetscObjectChangeTypeName((PetscObject)snes,SNESKSPTRANSPOSEONLY);CHKERRQ(ierr); 1391ef27442SStefano Zampini kspo = (SNES_KSPONLY*)snes->data; 1401ef27442SStefano Zampini kspo->transpose_solve = PETSC_TRUE; 141b79b07cfSJed Brown PetscFunctionReturn(0); 142b79b07cfSJed Brown } 143