OGS
anonymous_namespace{PETScNonlinearSolver.cpp} Namespace Reference

Classes

struct  PetscContext
 

Functions

PetscErrorCode updateResidual (SNES, Vec x, Vec petsc_r, void *petsc_context)
 
PetscErrorCode updateJacobian (SNES, Vec, Mat J, Mat, void *petsc_context)
 

Function Documentation

◆ updateJacobian()

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

Definition at line 63 of file PETScNonlinearSolver.cpp.

65 {
66  DBUG("PETScNonlinearSolver: jacobian callback called.");
67  // Assume the system is already assembled.
68  // Copy the results into petsc vectors.
69 
70  auto context = static_cast<PetscContext*>(petsc_context);
71  MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);
72  MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);
73  MatCopy(context->J->getRawMatrix(), J, DIFFERENT_NONZERO_PATTERN);
74  return 0;
75 }
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27

References DBUG().

Referenced by NumLib::PETScNonlinearSolver::solve().

◆ updateResidual()

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

Definition at line 33 of file PETScNonlinearSolver.cpp.

35 {
36  auto context = static_cast<PetscContext*>(petsc_context);
37 
38  DBUG("PETScNonlinearSolver: residual callback called.");
39 
40  VecCopy(x, context->x[context->process_id]->getRawVector());
41 
42  // Assemble in ogs context.
43  BaseLib::RunTime time_assembly;
44  time_assembly.start();
45  context->system->assemble(context->x, context->x_prev, context->process_id);
46 
47  INFO("[time] Assembly took {} s.", time_assembly.elapsed());
48  context->system->getResidual(*context->x[context->process_id],
49  *context->x_prev[context->process_id],
50  *context->r);
51  context->r->finalizeAssembly();
52  context->J->finalizeAssembly();
53 
54  context->system->getJacobian(*context->J);
55  context->system->applyKnownSolutionsNewton(
56  *context->J, *context->r, *context->x[context->process_id]);
57 
58  VecCopy(context->r->getRawVector(), petsc_r);
59 
60  return 0;
61 }
void INFO(char const *fmt, Args const &... args)
Definition: Logging.h:32
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().

Referenced by NumLib::PETScNonlinearSolver::solve().