FEI Package Browser (Single Doxygen Collection)
Version of the Day
Toggle main menu visibility
Loading...
Searching...
No Matches
base
snl_fei_LinearSystem_FEData.cpp
Go to the documentation of this file.
1
/*--------------------------------------------------------------------*/
2
/* Copyright 2005 Sandia Corporation. */
3
/* Under the terms of Contract DE-AC04-94AL85000, there is a */
4
/* non-exclusive license for use of this work by or on behalf */
5
/* of the U.S. Government. Export of this program may require */
6
/* a license from the United States Government. */
7
/*--------------------------------------------------------------------*/
8
#include <math.h>
9
10
#include <
fei_macros.hpp
>
11
12
#include <
fei_FiniteElementData.hpp
>
13
14
#include <
fei_CommUtils.hpp
>
15
#include <
snl_fei_LinearSystem_FEData.hpp
>
16
#include <
fei_VectorSpace.hpp
>
17
#include <
fei_MatrixGraph.hpp
>
18
#include <
fei_Matrix_Impl.hpp
>
19
#include <
snl_fei_Constraint.hpp
>
20
#include <
snl_fei_Utils.hpp
>
21
#include <
fei_impl_utils.hpp
>
22
23
#include <
fei_DirichletBCRecord.hpp
>
24
#include <
fei_DirichletBCManager.hpp
>
25
#include <
fei_EqnBuffer.hpp
>
26
#include <
fei_FEDataFilter.hpp
>
27
28
#undef fei_file
29
#define fei_file "snl_fei_LinearSystem_FEData.cpp"
30
#include <
fei_ErrMacros.hpp
>
31
32
//----------------------------------------------------------------------------
33
snl_fei::LinearSystem_FEData::LinearSystem_FEData
(
fei::SharedPtr<FiniteElementData>
& feData,
34
fei::SharedPtr<fei::MatrixGraph>
& matrixGraph)
35
:
fei
::
LinearSystem
(matrixGraph),
36
comm_
(matrixGraph->getRowSpace()->getCommunicator()),
37
localProc_
(0),
38
numProcs_
(1),
39
feData_
(feData),
40
lookup_
(NULL)
41
{
42
localProc_
=
fei::localProc
(
comm_
);
43
numProcs_
=
fei::numProcs
(
comm_
);
44
}
45
46
//----------------------------------------------------------------------------
47
snl_fei::LinearSystem_FEData::~LinearSystem_FEData
()
48
{
49
}
50
51
//----------------------------------------------------------------------------
52
bool
snl_fei::LinearSystem_FEData::eqnIsEssentialBC
(
int
globalEqnIndex)
const
53
{
54
fei::console_out
() <<
"LinearSystem_FEData::eqnIsEssentialBC NOT IMPLEMENTED!!"
<<
FEI_ENDL
;
55
return
(
false
);
56
}
57
58
//----------------------------------------------------------------------------
59
void
snl_fei::LinearSystem_FEData::getEssentialBCs
(std::vector<int>& bcEqns,
60
std::vector<double>& bcVals)
const
61
{
62
fei::console_out
() <<
"LinearSystem_FEData::getEssentialBC_Eqns NOT IMPLEMENTED!!"
<<
FEI_ENDL
;
63
}
64
65
//----------------------------------------------------------------------------
66
void
snl_fei::LinearSystem_FEData::getConstrainedEqns
(std::vector<int>& crEqns)
const
67
{
68
matrixGraph_
->getConstrainedIndices(crEqns);
69
}
70
71
//----------------------------------------------------------------------------
72
int
snl_fei::LinearSystem_FEData::loadComplete
(
bool
applyBCs,
73
bool
globalAssemble)
74
{
75
if
(
dbcManager_
== NULL) {
76
dbcManager_
=
new
fei::DirichletBCManager
(
matrixGraph_
->getRowSpace());
77
}
78
79
int
err = 0;
80
if
(
matrix_
.get() != NULL && globalAssemble) {
81
err =
matrix_
->gatherFromOverlap();
82
if
(err != 0) {
83
fei::console_out
() <<
"snl_fei::LinearSystem_FEData::loadComplete, ERROR in matrix."
84
<<
"gatherFromOverlap(), data may be incorrect."
<<
FEI_ENDL
;
85
}
86
}
87
88
fei::Barrier
(
comm_
);
89
90
if
(
rhs_
.get() != NULL && globalAssemble) {
91
err =
rhs_
->gatherFromOverlap();
92
if
(err != 0) {
93
fei::console_out
() <<
"snl_fei::LinearSystem_FEData::loadComplete, ERROR rhs."
94
<<
"gatherFromOverlap(), data may be incorrect."
<<
FEI_ENDL
;
95
}
96
}
97
98
fei::Barrier
(
comm_
);
99
100
CHK_ERR
(
implementBCs
(applyBCs) );
101
102
return
(
feData_
->loadComplete());
103
}
104
105
//----------------------------------------------------------------------------
106
int
snl_fei::LinearSystem_FEData::setBCValuesOnVector
(
fei::Vector
* vector)
107
{
108
ERReturn
(-1);
109
}
110
111
//----------------------------------------------------------------------------
112
int
snl_fei::LinearSystem_FEData::implementBCs
(
bool
applyBCs)
113
{
114
std::vector<int> essEqns;
115
std::vector<double> essGamma;
116
117
fei::SharedPtr<fei::FillableMat>
localBCeqns(
new
fei::FillableMat
);
118
fei::SharedPtr<fei::VectorSpace>
vecSpace =
matrixGraph_
->getRowSpace();
119
int
localsize = vecSpace->getNumIndices_Owned();
120
bool
zeroSharedRows =
false
;
121
fei::Matrix_Impl<fei::FillableMat>
bcEqns(localBCeqns,
matrixGraph_
, localsize, zeroSharedRows);
122
123
CHK_ERR
(
dbcManager_
->finalizeBCEqns(bcEqns) );
124
125
std::map<int,fei::FillableMat*>& remotes = bcEqns.
getRemotelyOwnedMatrices
();
126
std::map<int,fei::FillableMat*>::iterator
127
it = remotes.begin(),
128
it_end = remotes.end();
129
for
(; it!=it_end; ++it) {
130
fei::impl_utils::separate_BC_eqns
( *(it->second), essEqns, essGamma);
131
}
132
133
CHK_ERR
( bcEqns.
gatherFromOverlap
() );
134
135
fei::impl_utils::separate_BC_eqns
( *(bcEqns.
getMatrix
()), essEqns, essGamma);
136
137
int
len = essEqns.size();
138
essEqns.resize(len*3);
139
int
* nodeNumbers = &essEqns[0]+len;
140
int
* dof_ids = nodeNumbers+len;
141
142
if
(len > 0) {
143
int
* essEqnsPtr = &essEqns[0];
144
double
* gammaPtr = &essGamma[0];
145
146
for
(
int
i=0; i<len; ++i) {
147
int
eqn = essEqnsPtr[i];
148
nodeNumbers[i] =
lookup_
->getAssociatedNodeNumber(eqn);
149
int
fieldID =
lookup_
->getAssociatedFieldID(eqn);
150
int
base_eqn =
lookup_
->getEqnNumber(nodeNumbers[i], fieldID);
151
dof_ids[i]=vecSpace->getFieldDofMap().get_dof_id(fieldID, eqn-base_eqn);
152
}
153
154
if
(applyBCs) {
155
CHK_ERR
(
feData_
->setDirichletBCs(len, nodeNumbers, dof_ids, gammaPtr) );
156
}
157
}
158
159
return
(0);
160
}
161
162
//----------------------------------------------------------------------------
163
int
snl_fei::LinearSystem_FEData::
164
loadLagrangeConstraint
(
int
constraintID,
165
const
double
*weights,
166
double
rhsValue)
167
{
168
return
(-1);
169
}
170
171
//----------------------------------------------------------------------------
172
int
snl_fei::LinearSystem_FEData::
173
loadPenaltyConstraint
(
int
constraintID,
174
const
double
*weights,
175
double
penaltyValue,
176
double
rhsValue)
177
{
178
return
(-1);
179
}
fei::DirichletBCManager
Definition
fei_DirichletBCManager.hpp:25
fei::FillableMat
Definition
fei_FillableMat.hpp:20
fei::LinearSystem::matrixGraph_
fei::SharedPtr< fei::MatrixGraph > matrixGraph_
Definition
fei_LinearSystem.hpp:209
fei::LinearSystem::LinearSystem
LinearSystem(fei::SharedPtr< fei::MatrixGraph > &matrixGraph)
Definition
fei_LinearSystem.cpp:17
fei::LinearSystem::dbcManager_
fei::DirichletBCManager * dbcManager_
Definition
fei_LinearSystem.hpp:210
fei::Matrix_Impl
Definition
fei_Matrix_Impl.hpp:53
fei::Matrix_Impl::gatherFromOverlap
int gatherFromOverlap(bool accumulate=true)
Definition
fei_Matrix_Impl.hpp:945
fei::Matrix_Impl::getMatrix
fei::SharedPtr< T > getMatrix()
Definition
fei_Matrix_Impl.hpp:87
fei::Matrix_core::getRemotelyOwnedMatrices
std::map< int, FillableMat * > & getRemotelyOwnedMatrices()
Definition
fei_Matrix_core.cpp:87
fei::SharedPtr
Definition
fei_SharedPtr.hpp:65
fei::Vector
Definition
fei_Vector.hpp:57
snl_fei::LinearSystem_FEData::rhs_
fei::SharedPtr< fei::Vector > rhs_
Definition
snl_fei_LinearSystem_FEData.hpp:107
snl_fei::LinearSystem_FEData::getConstrainedEqns
void getConstrainedEqns(std::vector< int > &crEqns) const
Definition
snl_fei_LinearSystem_FEData.cpp:66
snl_fei::LinearSystem_FEData::loadComplete
int loadComplete(bool applyBCs=true, bool globalAssemble=true)
Definition
snl_fei_LinearSystem_FEData.cpp:72
snl_fei::LinearSystem_FEData::loadPenaltyConstraint
int loadPenaltyConstraint(int constraintID, const double *weights, double penaltyValue, double rhsValue)
Definition
snl_fei_LinearSystem_FEData.cpp:173
snl_fei::LinearSystem_FEData::comm_
MPI_Comm comm_
Definition
snl_fei_LinearSystem_FEData.hpp:102
snl_fei::LinearSystem_FEData::lookup_
Lookup * lookup_
Definition
snl_fei_LinearSystem_FEData.hpp:109
snl_fei::LinearSystem_FEData::setBCValuesOnVector
int setBCValuesOnVector(fei::Vector *vector)
Definition
snl_fei_LinearSystem_FEData.cpp:106
snl_fei::LinearSystem_FEData::~LinearSystem_FEData
virtual ~LinearSystem_FEData()
Definition
snl_fei_LinearSystem_FEData.cpp:47
snl_fei::LinearSystem_FEData::numProcs_
int numProcs_
Definition
snl_fei_LinearSystem_FEData.hpp:104
snl_fei::LinearSystem_FEData::eqnIsEssentialBC
bool eqnIsEssentialBC(int globalEqnIndex) const
Definition
snl_fei_LinearSystem_FEData.cpp:52
snl_fei::LinearSystem_FEData::LinearSystem_FEData
LinearSystem_FEData(fei::SharedPtr< FiniteElementData > &fedata, fei::SharedPtr< fei::MatrixGraph > &matrixGraph)
Definition
snl_fei_LinearSystem_FEData.cpp:33
snl_fei::LinearSystem_FEData::matrix_
fei::SharedPtr< fei::Matrix > matrix_
Definition
snl_fei_LinearSystem_FEData.hpp:105
snl_fei::LinearSystem_FEData::loadLagrangeConstraint
int loadLagrangeConstraint(int constraintID, const double *weights, double rhsValue)
Definition
snl_fei_LinearSystem_FEData.cpp:164
snl_fei::LinearSystem_FEData::localProc_
int localProc_
Definition
snl_fei_LinearSystem_FEData.hpp:103
snl_fei::LinearSystem_FEData::implementBCs
int implementBCs(bool applyBCs)
Definition
snl_fei_LinearSystem_FEData.cpp:112
snl_fei::LinearSystem_FEData::getEssentialBCs
void getEssentialBCs(std::vector< int > &bcEqns, std::vector< double > &bcVals) const
Definition
snl_fei_LinearSystem_FEData.cpp:59
snl_fei::LinearSystem_FEData::feData_
fei::SharedPtr< FiniteElementData > feData_
Definition
snl_fei_LinearSystem_FEData.hpp:108
fei_CommUtils.hpp
fei_DirichletBCManager.hpp
fei_DirichletBCRecord.hpp
fei_EqnBuffer.hpp
fei_ErrMacros.hpp
ERReturn
#define ERReturn(a)
Definition
fei_ErrMacros.hpp:37
CHK_ERR
#define CHK_ERR(a)
Definition
fei_ErrMacros.hpp:26
fei_FEDataFilter.hpp
fei_FiniteElementData.hpp
fei_MatrixGraph.hpp
fei_Matrix_Impl.hpp
fei_VectorSpace.hpp
fei_impl_utils.hpp
FEI_ENDL
#define FEI_ENDL
Definition
fei_iostream.hpp:34
fei_macros.hpp
fei::impl_utils::separate_BC_eqns
void separate_BC_eqns(const fei::FillableMat &mat, std::vector< int > &bcEqns, std::vector< double > &bcVals)
Definition
fei_impl_utils.cpp:258
fei
Definition
fei_ArrayUtils.hpp:16
fei::localProc
int localProc(MPI_Comm comm)
Definition
fei_CommUtils.cpp:13
fei::console_out
std::ostream & console_out()
Definition
fei_console_ostream.cpp:26
fei::numProcs
int numProcs(MPI_Comm comm)
Definition
fei_CommUtils.cpp:25
fei::Barrier
void Barrier(MPI_Comm comm)
Definition
fei_CommUtils.cpp:37
snl_fei_Constraint.hpp
snl_fei_LinearSystem_FEData.hpp
snl_fei_Utils.hpp
Generated by
1.17.0