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

Hotfix for Transient problems in parallel #70

Merged
merged 2 commits into from
Dec 19, 2024
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
4 changes: 3 additions & 1 deletion include/problem_operators/problem_operator_interface.h
Original file line number Diff line number Diff line change
Expand Up @@ -11,9 +11,11 @@ class ProblemOperatorInterface
virtual ~ProblemOperatorInterface() = default;

virtual void SetGridFunctions();
virtual void SetTestVariablesFromTrueVectors();
virtual void SetTrialVariablesFromTrueVectors();
virtual void Init(mfem::BlockVector & X);

mfem::Array<int> _true_offsets, _block_true_offsets;
mfem::Array<int> _block_true_offsets;

mfem::BlockVector _true_x, _true_rhs;
mfem::OperatorHandle _equation_system_operator;
Expand Down
3 changes: 3 additions & 0 deletions src/executioners/MFEMTransient.C
Original file line number Diff line number Diff line change
Expand Up @@ -52,6 +52,9 @@ MFEMTransient::step(double dt, int it) const
// Advance time step.
_problem_data._ode_solver->Step(_problem_data._f, _t, dt);

// Synchonise time dependent GridFunctions with updated DoF data.
_problem_operator->SetTestVariablesFromTrueVectors();

// Sync Host/Device
_problem_data._f.HostRead();

Expand Down
2 changes: 1 addition & 1 deletion src/problem_operators/problem_operator.C
Original file line number Diff line number Diff line change
Expand Up @@ -7,7 +7,7 @@ void
ProblemOperator::SetGridFunctions()
{
ProblemOperatorInterface::SetGridFunctions();
width = height = _true_offsets[_trial_variables.size()];
width = height = _block_true_offsets[_trial_variables.size()];
};

} // namespace platypus
30 changes: 20 additions & 10 deletions src/problem_operators/problem_operator_interface.C
Original file line number Diff line number Diff line change
Expand Up @@ -17,25 +17,35 @@ ProblemOperatorInterface::SetGridFunctions()
}
_block_true_offsets.PartialSum();

_true_offsets.SetSize(_trial_variables.size() + 1);
_true_offsets[0] = 0;
for (unsigned int ind = 0; ind < _trial_variables.size(); ++ind)
{
_true_offsets[ind + 1] = _trial_variables.at(ind)->ParFESpace()->GetVSize();
}
_true_offsets.PartialSum();

_true_x.Update(_block_true_offsets);
_true_rhs.Update(_block_true_offsets);
}

void
ProblemOperatorInterface::Init(mfem::BlockVector & X)
{
X.Update(_true_offsets);
X.Update(_block_true_offsets);
for (size_t i = 0; i < _test_variables.size(); ++i)
{
_test_variables.at(i)->MakeRef(_test_variables.at(i)->ParFESpace(), X, _true_offsets[i]);
_test_variables.at(i)->MakeTRef(_test_variables.at(i)->ParFESpace(), X, _block_true_offsets[i]);
}
}

void
ProblemOperatorInterface::SetTestVariablesFromTrueVectors()
{
for (unsigned int ind = 0; ind < _test_variables.size(); ++ind)
{
_test_variables.at(ind)->SetFromTrueVector();
}
}

void
ProblemOperatorInterface::SetTrialVariablesFromTrueVectors()
{
for (unsigned int ind = 0; ind < _trial_variables.size(); ++ind)
{
_trial_variables.at(ind)->SetFromTrueVector();
}
}

Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -27,10 +27,11 @@ TimeDomainEquationSystemProblemOperator::ImplicitSolve(const double dt,
mfem::Vector & dX_dt)
{
dX_dt = 0.0;
SetTestVariablesFromTrueVectors();
for (unsigned int ind = 0; ind < _trial_variables.size(); ++ind)
{
_trial_variables.at(ind)->MakeRef(
_trial_variables.at(ind)->ParFESpace(), dX_dt, _true_offsets[ind]);
_trial_variables.at(ind)->MakeTRef(
_trial_variables.at(ind)->ParFESpace(), dX_dt, _block_true_offsets[ind]);
}
const double time = GetTime();
for (auto & coef : _problem._scalar_manager)
Expand All @@ -50,6 +51,7 @@ TimeDomainEquationSystemProblemOperator::ImplicitSolve(const double dt,
_problem._nonlinear_solver->SetSolver(*_problem._jacobian_solver);
_problem._nonlinear_solver->SetOperator(*GetEquationSystem());
_problem._nonlinear_solver->Mult(_true_rhs, dX_dt);
SetTrialVariablesFromTrueVectors();
}

void
Expand Down
2 changes: 1 addition & 1 deletion src/problem_operators/time_domain_problem_operator.C
Original file line number Diff line number Diff line change
Expand Up @@ -17,7 +17,7 @@ void
TimeDomainProblemOperator::SetGridFunctions()
{
ProblemOperatorInterface::SetGridFunctions();
width = height = _true_offsets[_trial_variables.size()];
width = height = _block_true_offsets[_trial_variables.size()];
}

} // namespace platypus
Loading