v0.14.0
FormsIntegrators.cpp
Go to the documentation of this file.
1 /** \file FormsIntegrators.cpp
2 
3 \brief Implementation of Elements on Entities for Forces and Sources
4 
5 */
6 
7 
8 
9 namespace MoFEM {
10 
11 //! [Storage and set boundary conditions]
12 
13 
15 
16 OpSetBc::OpSetBc(std::string field_name, bool yes_set,
17  boost::shared_ptr<std::vector<unsigned char>> boundary_marker)
20  yesSet(yes_set), boundaryMarker(boundary_marker) {}
21 
22 MoFEMErrorCode OpSetBc::doWork(int side, EntityType type,
25  if (boundaryMarker) {
26  if (!data.getIndices().empty())
27  if (!data.getFieldEntities().empty()) {
28  if (auto e_ptr = data.getFieldEntities()[0]) {
29  auto indices = data.getIndices();
30  for (int r = 0; r != data.getIndices().size(); ++r) {
31  const auto loc_index = data.getLocalIndices()[r];
32  if (loc_index >= 0) {
33  if (yesSet) {
34  if ((*boundaryMarker)[loc_index]) {
35  indices[r] = -1;
36  }
37  } else {
38  if (!(*boundaryMarker)[loc_index]) {
39  indices[r] = -1;
40  }
41  }
42  }
43  }
44  EssentialBcStorage::feStorage[e_ptr->getName()].push_back(
45  boost::make_shared<EssentialBcStorage>(indices));
46  e_ptr->getWeakStoragePtr() =
47  EssentialBcStorage::feStorage[e_ptr->getName()].back();
48  }
49  }
50  }
52 }
53 
57 
58 MoFEMErrorCode OpUnSetBc::doWork(int side, EntityType type,
63 }
64 
65 /**
66  * @brief Set values to vector in operator
67  *
68  * @param V
69  * @param data
70  * @param ptr
71  * @param iora
72  * @return MoFEMErrorCode
73  */
74 template <>
77  const EntitiesFieldData::EntData &data,
78  const double *ptr, InsertMode iora) {
80 
81  if(!V)
82  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
83  "Pointer to PETSc vector is null");
84 
85  CHKERR VecSetOption(V, VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE);
86  if (!data.getFieldEntities().empty()) {
87  if (auto e_ptr = data.getFieldEntities()[0]) {
88  if (auto stored_data_ptr =
89  e_ptr->getSharedStoragePtr<EssentialBcStorage>()) {
90  return ::VecSetValues(V, stored_data_ptr->entityIndices.size(),
91  &*stored_data_ptr->entityIndices.begin(), ptr,
92  iora);
93  }
94  }
95  }
96  return ::VecSetValues(V, data.getIndices().size(),
97  &*data.getIndices().begin(), ptr, iora);
99 }
100 
101 /**
102  * @brief Set values to matrix in operator
103  *
104  * @param M
105  * @param row_data
106  * @param col_data
107  * @param ptr
108  * @param iora
109  * @return MoFEMErrorCode
110  */
111 template <>
113  Mat M, const EntitiesFieldData::EntData &row_data,
114  const EntitiesFieldData::EntData &col_data, const double *ptr,
115  InsertMode iora) {
116 
117  if(!M)
118  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
119  "Pointer to PETSc matrix is null");
120 
121  if (!row_data.getFieldEntities().empty()) {
122  if (auto e_ptr = row_data.getFieldEntities()[0]) {
123  if (auto stored_data_ptr =
124  e_ptr->getSharedStoragePtr<EssentialBcStorage>()) {
125  return ::MatSetValues(M, stored_data_ptr->entityIndices.size(),
126  &*stored_data_ptr->entityIndices.begin(),
127  col_data.getIndices().size(),
128  &*col_data.getIndices().begin(), ptr, iora);
129  }
130  }
131  }
133  M, row_data.getIndices().size(), &*row_data.getIndices().begin(),
134  col_data.getIndices().size(), &*col_data.getIndices().begin(), ptr, iora);
135 }
136 
137 //! [Storage and set boundary conditions
138 
139 
140 }
MoFEM::OpSetBc::OpSetBc
OpSetBc(std::string field_name, bool yes_set, boost::shared_ptr< std::vector< unsigned char >> boundary_marker)
Definition: FormsIntegrators.cpp:16
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
MoFEM::EntitiesFieldData::EntData::getFieldEntities
const VectorFieldEntities & getFieldEntities() const
get field entities
Definition: EntitiesFieldData.hpp:1277
MoFEM::MatSetValues
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
Definition: EntitiesFieldData.hpp:1644
MoFEM::EntitiesFieldData::EntData::getLocalIndices
const VectorInt & getLocalIndices() const
get local indices of dofs on entity
Definition: EntitiesFieldData.hpp:1229
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::EssentialBcStorage::HashVectorStorage
map< std::string, std::vector< boost::shared_ptr< EssentialBcStorage > >> HashVectorStorage
Store modifed indices by field.
Definition: FormsIntegrators.hpp:24
MoFEM::ForcesAndSourcesCore::UserDataOperator::rowFieldName
std::string rowFieldName
Definition: ForcesAndSourcesCore.hpp:577
MoFEM::OpUnSetBc::OpUnSetBc
OpUnSetBc(std::string field_name)
Definition: FormsIntegrators.cpp:54
MoFEM::VecSetValues
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
Definition: EntitiesFieldData.hpp:1589
sdf.r
int r
Definition: sdf.py:8
MoFEM::OpSetBc::boundaryMarker
boost::shared_ptr< std::vector< unsigned char > > boundaryMarker
Definition: FormsIntegrators.hpp:46
PlasticOps::M
FTensor::Index< 'M', 3 > M
Definition: PlasticOps.hpp:117
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::OpUnSetBc::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: FormsIntegrators.cpp:58
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::MatSetValues< EssentialBcStorage >
MoFEMErrorCode MatSetValues< EssentialBcStorage >(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Set values to matrix in operator.
Definition: FormsIntegrators.cpp:112
MoFEM::ForcesAndSourcesCore::UserDataOperator
Definition: ForcesAndSourcesCore.hpp:549
convert.type
type
Definition: convert.py:64
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1214
MoFEM::EssentialBcStorage::feStorage
static HashVectorStorage feStorage
[Storage and set boundary conditions]
Definition: FormsIntegrators.hpp:25
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MoFEM::OpSetBc::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: FormsIntegrators.cpp:22
MoFEM::ForcesAndSourcesCore
structure to get information form mofem into EntitiesFieldData
Definition: ForcesAndSourcesCore.hpp:22
MoFEM::EssentialBcStorage
[Storage and set boundary conditions]
Definition: FormsIntegrators.hpp:14
MoFEM::OpSetBc::yesSet
bool yesSet
Definition: FormsIntegrators.hpp:45
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::VecSetValues< EssentialBcStorage >
MoFEMErrorCode VecSetValues< EssentialBcStorage >(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Set values to vector in operator.
Definition: FormsIntegrators.cpp:76
MOFEM_INVALID_DATA
@ MOFEM_INVALID_DATA
Definition: definitions.h:36