Skip to content

Commit

Permalink
Refactor SteadyMeshAdaptiveSolver AMR output
Browse files Browse the repository at this point in the history
Allow for output after every AMR step via the 'output_every_amr' input
flag
  • Loading branch information
tradowsk committed Jan 17, 2018
1 parent e2cbb73 commit 641b0f1
Show file tree
Hide file tree
Showing 3 changed files with 81 additions and 51 deletions.
2 changes: 2 additions & 0 deletions src/solver/include/grins/steady_mesh_adaptive_solver.h
Original file line number Diff line number Diff line change
Expand Up @@ -60,6 +60,8 @@ namespace GRINS

protected:

bool amr_helper(SolverContext & context, unsigned int r_step, bool do_amr);

virtual void init_time_solver( MultiphysicsSystem* system );

void check_qoi_error_option_consistency( SolverContext& context );
Expand Down
125 changes: 74 additions & 51 deletions src/solver/src/steady_mesh_adaptive_solver.C
Original file line number Diff line number Diff line change
Expand Up @@ -85,6 +85,8 @@ namespace GRINS
libMesh::MeshBase& mesh = context.equation_system->get_mesh();
this->build_mesh_refinement( mesh );

bool converged = false;

/*! \todo This output cannot be toggled in the input file, but it should be able to be. */
std::cout << "==========================================================" << std::endl
<< "Performing " << _mesh_adaptivity_options.max_refinement_steps()
Expand All @@ -93,71 +95,92 @@ namespace GRINS

for ( unsigned int r_step = 0; r_step < _mesh_adaptivity_options.max_refinement_steps(); r_step++ )
{
std::cout << "==========================================================" << std::endl
std::cout << std::endl
<< "==========================================================" << std::endl
<< "Adaptive Refinement Step " << r_step << std::endl
<< "==========================================================" << std::endl;

// Solve the forward problem
context.system->solve();
converged = this->amr_helper(context,r_step,true);
if (converged)
break;

if( context.output_vis )
{
context.postprocessing->update_quantities( *(context.equation_system) );
context.vis->output( context.equation_system );
}
} // r_step for-loop

// Solve adjoint system
if(context.do_adjoint_solve)
this->steady_adjoint_solve(context);
// this will print output the error estimates and/or meshes after the final AMR run
if (! converged)
{
std::cout << std::endl
<< "==========================================================" << std::endl
<< "Post-AMR Output " << std::endl
<< "==========================================================" << std::endl;
this->amr_helper(context,_mesh_adaptivity_options.max_refinement_steps(),false);
}

if(context.output_adjoint)
context.vis->output_adjoint(context.equation_system, context.system);
return;
}

if( context.output_residual )
{
context.vis->output_residual( context.equation_system, context.system );
}
bool SteadyMeshAdaptiveSolver::amr_helper(SolverContext & context, unsigned int r_step, bool do_amr)
{
// Solve the forward problem
context.system->solve();

// Now we construct the data structures for the mesh refinement process
libMesh::ErrorVector error;
this->estimate_error_for_amr( context, error );
if( context.output_vis )
{
context.postprocessing->update_quantities( *(context.equation_system) );

// Get the global error estimate if you can and are asked to
if( _error_estimator_options.compute_qoi_error_estimate() )
for(unsigned int i = 0; i != context.system->qoi.size(); i++)
{
libMesh::AdjointRefinementEstimator* adjoint_ref_error_estimator = libMesh::cast_ptr<libMesh::AdjointRefinementEstimator*>( context.error_estimator.get() );
std::cout<<"The error estimate for QoI("<<i<<") is: "<<adjoint_ref_error_estimator->get_global_QoI_error_estimate(i)<<std::endl;
}
if (context.output_every_amr)
context.vis->output_amr( context.equation_system,r_step );
else
context.vis->output( context.equation_system );
}

// Check for convergence of error
bool converged = this->check_for_convergence( context, error );
// Solve adjoint system
if(context.do_adjoint_solve)
this->steady_adjoint_solve(context);

if( converged )
{
// Break out of adaptive loop
std::cout << "==========================================================" << std::endl
<< "Convergence detected!" << std::endl
<< "==========================================================" << std::endl;
break;
}
else
if(context.output_adjoint)
context.vis->output_adjoint(context.equation_system, context.system);

if( context.output_residual )
{
context.vis->output_residual( context.equation_system, context.system );
}

// Now we construct the data structures for the mesh refinement process
libMesh::ErrorVector error;
this->estimate_error_for_amr( context, error );

// Get the global error estimate if you can and are asked to
if( _error_estimator_options.compute_qoi_error_estimate() )
for(unsigned int i = 0; i != context.system->qoi.size(); i++)
{
libMesh::AdjointRefinementEstimator* adjoint_ref_error_estimator = libMesh::cast_ptr<libMesh::AdjointRefinementEstimator*>( context.error_estimator.get() );
std::cout<<"The error estimate for QoI("<<i<<") is: "<<adjoint_ref_error_estimator->get_global_QoI_error_estimate(i)<<std::endl;
}

// Check for convergence of error
bool converged = this->check_for_convergence(context,error);

if( converged )
{
// Break out of adaptive loop
std::cout << "==========================================================" << std::endl
<< "Convergence detected!" << std::endl
<< "==========================================================" << std::endl;
}
else
{
if (do_amr)
{
// Only bother refining if we're on the last step.
if( r_step < _mesh_adaptivity_options.max_refinement_steps() -1 )
{
this->perform_amr(context, error);

// It's helpful to print the qoi along the way, but only do it if the user
// asks for it
if( context.qoi_output->output_qoi_set() )
this->print_qoi(context);
}
// print the QoI before AMR
if( context.qoi_output->output_qoi_set() )
this->print_qoi(context,context.output_every_amr);

this->perform_amr(context, error);
}
}

} // r_step for-loop

return;
return converged;
}

void SteadyMeshAdaptiveSolver::adjoint_qoi_parameter_sensitivity
Expand Down
5 changes: 5 additions & 0 deletions src/solver/src/unsteady_mesh_adaptive_solver.C
Original file line number Diff line number Diff line change
Expand Up @@ -44,6 +44,11 @@ namespace GRINS

void UnsteadyMeshAdaptiveSolver::solve( SolverContext& context )
{
if (context.output_every_amr)
std::cout << std::endl
<< "WARNING: 'output_every_amr' flag is not yet implemented"
<< "in UnsteadyMeshAdaptiveSolver" <<std::endl;

context.system->deltat = this->_deltat;

libMesh::Real sim_time;
Expand Down

0 comments on commit 641b0f1

Please sign in to comment.