NeoN
WIP Prototype of a modern OpenFOAM core
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"
12
14
15// TODO decouple from fvcc
17
18
20{
21
22template<typename SolutionVectorType>
24 public TimeIntegratorBase<SolutionVectorType>::template Register<
25 BackwardEuler<SolutionVectorType>>
26{
27
28public:
29
30 using ValueType = typename SolutionVectorType::VectorValueType;
33
34 BackwardEuler(const Dictionary& schemeDict, const Dictionary& solutionDict)
35 : Base(schemeDict, solutionDict)
36 {}
37
38 static std::string name() { return "backwardEuler"; }
39
40 static std::string doc() { return "first order time integration method"; }
41
42 static std::string schema() { return "none"; }
43
44 void solve(
46 SolutionVectorType& solutionVector,
47 [[maybe_unused]] scalar t,
48 scalar dt
49 ) override
50 {
51 auto source = eqn.explicitOperation(solutionVector.size());
52 SolutionVectorType& oldSolutionVector = finiteVolume::cellCentred::oldTime(solutionVector);
53
54 // solutionVector.internalVector() = oldSolutionVector.internalVector() - source * dt;
55 // solutionVector.correctBoundaryConditions();
56 // solve sparse matrix system
57 // using ValueType = typename SolutionVectorType::ElementType;
58
59 // TODO decouple from fvcc specific implementation
60 auto sparsity = NeoN::finiteVolume::cellCentred::SparsityPattern(solutionVector.mesh());
65
66 eqn.implicitOperation(ls);
67
68 auto values = ls.matrix().values();
69 eqn.implicitOperation(ls, t, dt);
70
71 auto solver = NeoN::la::Solver(solutionVector.exec(), this->solutionDict_);
72 solver.solve(ls, solutionVector.internalVector());
73
74 // check if executor is GPU
75 if (std::holds_alternative<NeoN::GPUExecutor>(eqn.exec()))
76 {
77 Kokkos::fence();
78 }
79 oldSolutionVector.internalVector() = solutionVector.internalVector();
80 };
81
82 std::unique_ptr<TimeIntegratorBase<SolutionVectorType>> clone() const override
83 {
84 return std::make_unique<BackwardEuler>(*this);
85 }
86};
87
88
89} // 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.
VectorType & oldTime(VectorType &field)
Retrieves the old time field of a given field.
LinearSystem< ValueType, IndexType > createEmptyLinearSystem(const SparsityType &sparsity)
int32_t localIdx
Definition label.hpp:30
float scalar
Definition scalar.hpp:14