OGS
anonymous_namespace{PETScNonlinearSolver.cpp} Namespace Reference

Classes

struct  PetscContext

Functions

void printConvergenceReason (SNES const &snes_solver, SNESConvergedReason const &reason)
PetscErrorCode updateResidual (SNES, Vec x, Vec petsc_r, void *petsc_context)
PetscErrorCode updateJacobian (SNES, Vec, Mat J, Mat, void *petsc_context)

Function Documentation

◆ printConvergenceReason()

void anonymous_namespace{PETScNonlinearSolver.cpp}::printConvergenceReason ( SNES const & snes_solver,
SNESConvergedReason const & reason )

Definition at line 26 of file PETScNonlinearSolver.cpp.

28{
29 const char* strreason;
30 SNESGetConvergedReasonString(snes_solver, &strreason);
31 if (reason < 0)
32 {
33 INFO("PETScSNES diverged reason: {:s}.", strreason);
34 }
35 else
36 {
37 INFO("PETScSNES converged reason: {:s}.", strreason);
38 }
39}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28

References INFO().

◆ updateJacobian()

PetscErrorCode anonymous_namespace{PETScNonlinearSolver.cpp}::updateJacobian ( SNES ,
Vec ,
Mat J,
Mat ,
void * petsc_context )

Definition at line 71 of file PETScNonlinearSolver.cpp.

73{
74 DBUG("PETScNonlinearSolver: jacobian callback called.");
75 // Assume the system is already assembled.
76 // Copy the results into petsc vectors.
77
78 auto context = static_cast<PetscContext*>(petsc_context);
79 MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);
80 MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);
81 MatCopy(context->J->getRawMatrix(), J, DIFFERENT_NONZERO_PATTERN);
82 return 0;
83}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22

References DBUG().

◆ updateResidual()

PetscErrorCode anonymous_namespace{PETScNonlinearSolver.cpp}::updateResidual ( SNES ,
Vec x,
Vec petsc_r,
void * petsc_context )

Definition at line 41 of file PETScNonlinearSolver.cpp.

43{
44 auto context = static_cast<PetscContext*>(petsc_context);
45
46 DBUG("PETScNonlinearSolver: residual callback called.");
47
48 VecCopy(x, context->x[context->process_id]->getRawVector());
49
50 // Assemble in ogs context.
51 BaseLib::RunTime time_assembly;
52 time_assembly.start();
53 context->system->assemble(context->x, context->x_prev, context->process_id);
54
55 INFO("[time] Assembly took {} s.", time_assembly.elapsed());
56 context->system->getResidual(*context->x[context->process_id],
57 *context->x_prev[context->process_id],
58 *context->r);
59 context->r->finalizeAssembly();
60 context->J->finalizeAssembly();
61
62 context->system->getJacobian(*context->J);
63 context->system->applyKnownSolutionsPETScSNES(
64 *context->J, *context->r, *context->x[context->process_id]);
65
66 VecCopy(context->r->getRawVector(), petsc_r);
67
68 return 0;
69}
Count the running time.
Definition RunTime.h:18
double elapsed() const
Get the elapsed time in seconds.
Definition RunTime.h:31
void start()
Start the timer.
Definition RunTime.h:21

References DBUG(), BaseLib::RunTime::elapsed(), INFO(), and BaseLib::RunTime::start().