Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

PETSc stages #8

Merged
merged 8 commits into from
Dec 17, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion Project.toml
Original file line number Diff line number Diff line change
@@ -1,7 +1,7 @@
name = "LaMEM_C"
uuid = "29b6043d-3311-4c16-b469-b02bf85cb680"
authors = ["Anton Popov <[email protected]>", "Boris Kaus <[email protected]>"]
version = "2.1.1"
version = "2.1.2"

[deps]
CairoMakie = "13f3f980-e62b-5c42-98c6-ff1f3baf88f0"
Expand Down
3 changes: 1 addition & 2 deletions src/LaMEM.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -22,7 +22,6 @@ int main(int argc, char **argv)

// Initialize PETSC
ierr = PetscInitialize(&argc,&argv,(char *)0, help); CHKERRQ(ierr);

ModParam IOparam;
char str[_str_len_];

Expand All @@ -45,7 +44,7 @@ int main(int argc, char **argv)
if(IOparam.use == 0)
{
// Forward simulation
ierr = LaMEMLibMain(NULL); CHKERRQ(ierr);
ierr = LaMEMLibMain(NULL,IOparam.stages); CHKERRQ(ierr);
}
else
{
Expand Down
2 changes: 1 addition & 1 deletion src/LaMEM.h
Original file line number Diff line number Diff line change
Expand Up @@ -189,7 +189,7 @@ typedef pair <PetscInt, PetscInt> ipair;

// LaMEM library main function

PetscErrorCode LaMEMLibMain(void *param);
PetscErrorCode LaMEMLibMain(void *param,PetscLogStage stages[4]);

//-----------------------------------------------------------------------------
#endif
25 changes: 20 additions & 5 deletions src/LaMEMLib.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -42,7 +42,7 @@
#include "passive_tracer.h"

//---------------------------------------------------------------------------
PetscErrorCode LaMEMLibMain(void *param)
PetscErrorCode LaMEMLibMain(void *param,PetscLogStage stages[4])
{
LaMEMLib lm;
RunMode mode;
Expand All @@ -60,7 +60,7 @@ PetscErrorCode LaMEMLibMain(void *param)
PetscPrintf(PETSC_COMM_WORLD,"-------------------------------------------------------------------------- \n");
PetscPrintf(PETSC_COMM_WORLD," Lithosphere and Mantle Evolution Model \n");
PetscPrintf(PETSC_COMM_WORLD," Compiled: Date: %s - Time: %s \n",__DATE__,__TIME__ );
PetscPrintf(PETSC_COMM_WORLD," Version : 2.1.0 \n");
PetscPrintf(PETSC_COMM_WORLD," Version : 2.1.2 \n");
PetscPrintf(PETSC_COMM_WORLD,"-------------------------------------------------------------------------- \n");
PetscPrintf(PETSC_COMM_WORLD," STAGGERED-GRID FINITE DIFFERENCE CANONICAL IMPLEMENTATION \n");
PetscPrintf(PETSC_COMM_WORLD,"-------------------------------------------------------------------------- \n");
Expand Down Expand Up @@ -130,7 +130,7 @@ PetscErrorCode LaMEMLibMain(void *param)
else if(mode == _NORMAL_ || mode == _RESTART_)
{
// solve coupled nonlinear equations
ierr = LaMEMLibSolve(&lm, param); CHKERRQ(ierr);
ierr = LaMEMLibSolve(&lm, param,stages); CHKERRQ(ierr);
}

// destroy library objects
Expand Down Expand Up @@ -604,7 +604,7 @@ PetscErrorCode LaMEMLibSaveOutput(LaMEMLib *lm)
PetscFunctionReturn(0);
}
//---------------------------------------------------------------------------
PetscErrorCode LaMEMLibSolve(LaMEMLib *lm, void *param)
PetscErrorCode LaMEMLibSolve(LaMEMLib *lm, void *param, PetscLogStage stages[4])
{
PMat pm; // preconditioner matrix (to be removed!)
PCStokes pc; // Stokes preconditioner (to be removed!)
Expand All @@ -614,6 +614,7 @@ PetscErrorCode LaMEMLibSolve(LaMEMLib *lm, void *param)
PetscInt restart;
PetscLogDouble t;


PetscErrorCode ierr;
PetscFunctionBeginUser;

Expand Down Expand Up @@ -649,17 +650,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 +693,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 +716,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 +736,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
2 changes: 1 addition & 1 deletion src/LaMEMLib.h
Original file line number Diff line number Diff line change
Expand Up @@ -68,7 +68,7 @@ PetscErrorCode LaMEMLibSetLinks(LaMEMLib *lm);

PetscErrorCode LaMEMLibSaveOutput(LaMEMLib *lm, PetscInt dirInd);

PetscErrorCode LaMEMLibSolve(LaMEMLib *lm, void *param);
PetscErrorCode LaMEMLibSolve(LaMEMLib *lm, void *param, PetscLogStage stages[4]);

PetscErrorCode LaMEMLibDryRun(LaMEMLib *lm);

Expand Down
10 changes: 5 additions & 5 deletions src/adjoint.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -984,7 +984,7 @@ PetscErrorCode LaMEMAdjointMain(ModParam *IOparam)
// Only compute a forward model and the corresponding misfit

IOparam->BruteForce_FD = PETSC_TRUE; // to return early w/out cmomputing FD gradients
ierr = LaMEMLibMain(IOparam); CHKERRQ(ierr);
ierr = LaMEMLibMain(IOparam,IOparam->stages); CHKERRQ(ierr);

// Print overview of cost function & gradients
ierr = PrintCostFunction(IOparam); CHKERRQ(ierr);
Expand Down Expand Up @@ -1079,7 +1079,7 @@ PetscErrorCode LaMEMAdjointMain(ModParam *IOparam)
else if(IOparam->use == _syntheticforwardrun_)
{
// call LaMEM main library function
ierr = LaMEMLibMain(IOparam); CHKERRQ(ierr);
ierr = LaMEMLibMain(IOparam,IOparam->stages); CHKERRQ(ierr);

// Save output
PetscViewer viewerVel;
Expand Down Expand Up @@ -1212,7 +1212,7 @@ PetscErrorCode ComputeGradientsAndObjectiveFunction(Vec Parameters, PetscScalar

// Adjoint gradients: Call LaMEM main library function once (computes gradients @ the end)
IOparam->BruteForce_FD = PETSC_FALSE;
ierr = LaMEMLibMain(IOparam); CHKERRQ(ierr);
ierr = LaMEMLibMain(IOparam,IOparam->stages); CHKERRQ(ierr);

// Print overview of cost function & gradients
ierr = PrintCostFunction(IOparam); CHKERRQ(ierr);
Expand Down Expand Up @@ -1827,7 +1827,7 @@ PetscErrorCode AdjointFiniteDifferenceGradients(ModParam *IOparam)
if (FD_Adjoint){

// Call LaMEM
ierr = LaMEMLibMain(IOparam); CHKERRQ(ierr); // call LaMEM
ierr = LaMEMLibMain(IOparam,IOparam->stages); CHKERRQ(ierr); // call LaMEM
Misfit_ref = IOparam->mfit;

PetscPrintf(PETSC_COMM_WORLD,"| ************************************************************************ \n");
Expand Down Expand Up @@ -1863,7 +1863,7 @@ PetscErrorCode AdjointFiniteDifferenceGradients(ModParam *IOparam)
ierr = CopyParameterToLaMEMCommandLine(IOparam, CurVal + Perturb, j); CHKERRQ(ierr);

// Compute solution with updated parameter
ierr = LaMEMLibMain(IOparam); CHKERRQ(ierr);
ierr = LaMEMLibMain(IOparam,IOparam->stages); CHKERRQ(ierr);
Misfit_pert = IOparam->mfit;

// FD gradient
Expand Down
2 changes: 2 additions & 0 deletions src/objFunct.h
Original file line number Diff line number Diff line change
Expand Up @@ -115,6 +115,8 @@ struct ModParam
PetscScalar Avel_num[_MAX_OBS_]; // Numerically computed velocity at the comparison points
PetscBool Apoint_on_proc[_MAX_OBS_]; // Is the observation point on the current processor or not (simplified printing)?
char ScalLawFilename[_str_len_]; // Name of scaling law file

PetscLogStage stages[4]; /* Create stages for PETSC profiling */
};

// observation type
Expand Down