NeoN
A framework for CFD software
Loading...
Searching...
No Matches
backwardEuler.hpp
Go to the documentation of this file.
1// SPDX-License-Identifier: MIT
2//
3// SPDX-FileCopyrightText: 2023 NeoN authors
4
5#pragma once
6
11#include "NeoN/dsl/solver.hpp"
14
15
17{
18
19template<typename SolutionVectorType>
21 public TimeIntegratorBase<SolutionVectorType>::template Register<
22 BackwardEuler<SolutionVectorType>>
23{
24
25public:
26
27 using ValueType = typename SolutionVectorType::VectorValueType;
30
31 BackwardEuler(const Dictionary& schemeDict, const Dictionary& solutionDict)
32 : Base(schemeDict, solutionDict)
33 {}
34
35 static std::string name() { return "backwardEuler"; }
36
37 static std::string doc() { return "first order time integration method"; }
38
39 static std::string schema() { return "none"; }
40
41 void solve(
42 dsl::Expression<ValueType>& eqn, SolutionVectorType& solutionVector, scalar t, scalar dt
43 ) override
44 {
45 auto source = eqn.explicitOperation(solutionVector.size());
46 auto sparsity = la::SparsityPattern(solutionVector.mesh());
47 auto ls = la::createEmptyLinearSystem<ValueType, localIdx>(solutionVector.mesh(), sparsity);
48
49 eqn.implicitOperation(ls); // add spatial operators
50 eqn.implicitOperation(ls, t, dt); // add temporal operators
51
52 auto solver = NeoN::la::Solver(solutionVector.exec(), this->solutionDict_);
53 solver.solve(ls, solutionVector.internalVector());
54 // check if executor is GPU
55 if (std::holds_alternative<NeoN::GPUExecutor>(eqn.exec()))
56 {
57 Kokkos::fence();
58 }
59 };
60
61 std::unique_ptr<TimeIntegratorBase<SolutionVectorType>> clone() const override
62 {
63 return std::make_unique<BackwardEuler>(*this);
64 }
65};
66
67
68} // namespace NeoN
A class representing a dictionary that stores key-value pairs.
void implicitOperation(la::LinearSystem< ValueType, localIdx > &ls)
Vector< ValueType > explicitOperation(localIdx nCells) const
const Executor & exec() const
BackwardEuler(const Dictionary &schemeDict, const Dictionary &solutionDict)
TimeIntegratorBase< SolutionVectorType >::template Register< BackwardEuler< SolutionVectorType > > Base
std::unique_ptr< TimeIntegratorBase< SolutionVectorType > > clone() const override
void solve(dsl::Expression< ValueType > &eqn, SolutionVectorType &solutionVector, scalar t, scalar dt) override
typename SolutionVectorType::VectorValueType ValueType
A template class for registering derived classes with a base class.
float scalar
Definition scalar.hpp:14