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 33 of file PETScNonlinearSolver.cpp.

35{
36 const char* strreason;
37 SNESGetConvergedReasonString(snes_solver, &strreason);
38 if (reason < 0)
39 {
40 INFO("PETScSNES diverged reason: {:s}.", strreason);
41 }
42 else
43 {
44 INFO("PETScSNES converged reason: {:s}.", strreason);
45 }
46}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35

References INFO().

◆ updateJacobian()

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

Definition at line 78 of file PETScNonlinearSolver.cpp.

80{
81 DBUG("PETScNonlinearSolver: jacobian callback called.");
82 // Assume the system is already assembled.
83 // Copy the results into petsc vectors.
84
85 auto context = static_cast<PetscContext*>(petsc_context);
86 MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);
87 MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);
88 MatCopy(context->J->getRawMatrix(), J, DIFFERENT_NONZERO_PATTERN);
89 return 0;
90}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30

References DBUG().

◆ updateResidual()

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

Definition at line 48 of file PETScNonlinearSolver.cpp.

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

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