Skip to content

Commit

Permalink
Added 4 stages
Browse files Browse the repository at this point in the history
  • Loading branch information
IskanderI committed Dec 10, 2023
1 parent 7694caa commit 2c074ea
Showing 1 changed file with 23 additions and 1 deletion.
24 changes: 23 additions & 1 deletion src/LaMEMLib.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -613,6 +613,8 @@ PetscErrorCode LaMEMLibSolve(LaMEMLib *lm, void *param)
SNES snes; // PETSc nonlinear solver
PetscInt restart;
PetscLogDouble t;
PetscLogStage stages[4]; /* Create stages for PETSC profiling */


PetscErrorCode ierr;
PetscFunctionBeginUser;
Expand All @@ -622,6 +624,12 @@ PetscErrorCode LaMEMLibSolve(LaMEMLib *lm, void *param)
ierr = PCStokesCreate(&pc, pm); CHKERRQ(ierr);
ierr = NLSolCreate(&nl, pc, &snes); CHKERRQ(ierr);

/* Register names of the stages*/
PetscCall(PetscLogStageRegister("Thermal solver", &stages[0]));
PetscCall(PetscLogStageRegister("SNES solve", &stages[1]));
PetscCall(PetscLogStageRegister("Advect markers", &stages[2]));
PetscCall(PetscLogStageRegister("I/O", &stages[3]));

//==============
// INITIAL GUESS
//==============
Expand Down Expand Up @@ -649,17 +657,23 @@ PetscErrorCode LaMEMLibSolve(LaMEMLib *lm, void *param)
// initialize boundary constraint vectors
ierr = BCApply(&lm->bc); CHKERRQ(ierr);

PetscCall(PetscLogStagePush(stages[0])); /* Start profiling stage*/

// initialize temperature
ierr = JacResInitTemp(&lm->jr); CHKERRQ(ierr);

PetscCall(PetscLogStagePop()); /* Stop profiling stage*/
// compute elastic parameters
ierr = JacResGetI2Gdt(&lm->jr); CHKERRQ(ierr);

// solve nonlinear equation system with SNES
PetscTime(&t);

PetscCall(PetscLogStagePush(stages[1])); /* Start profiling stage*/

ierr = SNESSolve(snes, NULL, lm->jr.gsol); CHKERRQ(ierr);

PetscCall(PetscLogStagePop()); /* Stop profiling stage*/
// print analyze convergence/divergence reason & iteration count
ierr = SNESPrintConvergedReason(snes, t); CHKERRQ(ierr);

Expand All @@ -686,6 +700,8 @@ PetscErrorCode LaMEMLibSolve(LaMEMLib *lm, void *param)
// MARKER & FREE SURFACE ADVECTION + EROSION
//==========================================

PetscCall(PetscLogStagePush(stages[2])); /* Start profiling stage*/

// calculate current time step
ierr = ADVSelectTimeStep(&lm->actx, &restart); CHKERRQ(ierr);

Expand All @@ -707,6 +723,8 @@ PetscErrorCode LaMEMLibSolve(LaMEMLib *lm, void *param)
// Advect Passive tracers
ierr = ADVAdvectPassiveTracer(&lm->actx); CHKERRQ(ierr);

PetscCall(PetscLogStagePop()); /* Stop profiling stage*/

// apply erosion to the free surface
ierr = FreeSurfAppErosion(&lm->surf); CHKERRQ(ierr);

Expand All @@ -725,10 +743,14 @@ PetscErrorCode LaMEMLibSolve(LaMEMLib *lm, void *param)

// update time stamp and counter
ierr = TSSolStepForward(&lm->ts); CHKERRQ(ierr);


PetscCall(PetscLogStagePush(stages[3])); /* Start profiling stage*/

// grid & marker output
ierr = LaMEMLibSaveOutput(lm); CHKERRQ(ierr);

PetscCall(PetscLogStagePop()); /* Stop profiling stage*/

// restart database
ierr = LaMEMLibSaveRestart(lm); CHKERRQ(ierr);

Expand Down

0 comments on commit 2c074ea

Please sign in to comment.