-
Notifications
You must be signed in to change notification settings - Fork 8
/
Copy pathgsBaseAssembler.h
111 lines (83 loc) · 4.54 KB
/
gsBaseAssembler.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
/** @file gsBaseAssembler.h
@brief Base class for assemblers of gsElasticity.
This file is part of the G+Smo library.
This Source Code Form is subject to the terms of the Mozilla Public
License, v. 2.0. If a copy of the MPL was not distributed with this
file, You can obtain one at http://mozilla.org/MPL/2.0/.
Author(s):
A.Shamanskiy (2016 - ...., TU Kaiserslautern)
*/
#pragma once
#include <gsAssembler/gsAssembler.h>
namespace gismo
{
/** @brief Extends the gsAssembler class by adding functionality necessary for a general nonlinear solver.
* Potentially, can be merged back into gsAssembler.
*/
template <class T>
class gsBaseAssembler : public gsAssembler<T>
{
public:
typedef memory::shared_ptr<gsBaseAssembler> Ptr;
typedef memory::unique_ptr<gsBaseAssembler> uPtr;
/// Assembles the tangential linear system for Newton's method given the current solution
/// in the form of free and fixed/Dirichelt degrees of freedom.
/// Checks if the current solution is valid (Newton's solver can exit safely if invalid).
virtual bool assemble(const gsMatrix<T> & solutionVector,
const std::vector<gsMatrix<T> > & fixedDDoFs) = 0;
/// assembly procedure for linear problems
virtual void assemble(bool /* saveEliminationMatrix */) {};
virtual void assemble() { assemble(false); };
virtual void assemble(const gsMultiPatch<T> & /* curSolution */)
{GISMO_NO_IMPLEMENTATION}
/// Returns number of free degrees of freedom
virtual int numDofs() const { return gsAssembler<T>::numDofs(); }
/// Constructs solution as a gsMultiPatch object from the solution vector and fixed DoFs
virtual void constructSolution(const gsMatrix<T> & solVector,
const std::vector<gsMatrix<T> > & fixedDDofs,
gsMultiPatch<T> & result,
const gsVector<index_t> & unknowns) const;
virtual void constructSolution(const gsMatrix<T>& /* solVector */,
gsMultiPatch<T>& /* result */, short_t /* unk */ = 0) const
{GISMO_NO_IMPLEMENTATION}
virtual void constructSolution(const gsMatrix<T>& /* solVector */,
gsMultiPatch<T>& /* result */,
const gsVector<index_t> & /* unknowns */) const
{GISMO_NO_IMPLEMENTATION}
virtual void constructSolution(const gsMatrix<T> & /* solVector */,
const std::vector<gsMatrix<T> > & /* fixedDDofs */,
gsMultiPatch<T> & /* result */) const {};
//--------------------- DIRICHLET BC SHENANIGANS ----------------------------------//
/** @brief Set Dirichet degrees of freedom on a given side of a given patch from a given matrix.
*
* A usual source of degrees of freedom is another geometry where it is known that the corresponding
* bases match. The function is not safe in that it assumes a correct numbering of DoFs in a given
* matrix. To allocate space for these DoFs in the assembler, add an empty/zero Dirichlet boundary condition
* to gsBoundaryCondtions container that is passed to the assembler constructor.
*/
virtual void setFixedDofs(index_t patch, boxSide side, const gsMatrix<T> & ddofs, bool oneUnk = false);
/// set all fixed degrees of freedom
virtual void setFixedDofs(const std::vector<gsMatrix<T> > & ddofs);
/// get fixed degrees of freedom corresponding to a given part of the bdry.
/// each column of the resulting matrix correspond to one variable/component of the vector-valued vairable
virtual void getFixedDofs(index_t patch, boxSide side, gsMatrix<T> & ddofs) const;
/// get the size of the Dirichlet vector for elimination
virtual index_t numFixedDofs() const;
/// @brief Eliminates new Dirichelt degrees of fredom
virtual void eliminateFixedDofs();
//virtual void modifyDirichletDofs(size_t patch, boxSide side, const gsMatrix<T> & ddofs);
//--------------------- OTHER ----------------------------------//
virtual void setRHS(const gsMatrix<T> & rhs) {m_system.rhs() = rhs;}
virtual void setMatrix(const gsSparseMatrix<T> & matrix) {m_system.matrix() = matrix;}
protected:
using gsAssembler<T>::m_pde_ptr;
using gsAssembler<T>::m_bases;
using gsAssembler<T>::m_system;
using gsAssembler<T>::m_ddof;
gsSparseMatrix<T> eliminationMatrix;
gsMatrix<T> rhsWithZeroDDofs;
};
} // namespace ends
#ifndef GISMO_BUILD_LIB
#include GISMO_HPP_HEADER(gsBaseAssembler.hpp)
#endif