v0.13.2
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
OpRotatingFrameTools::OpVolumeSideAssembleRowData Struct Reference

#include <users_modules/multifield_plasticity/src/RotatingFrameOperators.hpp>

Inheritance diagram for OpRotatingFrameTools::OpVolumeSideAssembleRowData:
[legend]
Collaboration diagram for OpRotatingFrameTools::OpVolumeSideAssembleRowData:
[legend]

Public Member Functions

 OpVolumeSideAssembleRowData (const std::string field_name, boost::shared_ptr< CommonData > common_data_ptr)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 
MoFEMErrorCode iNtegrate (EntData &row_data, EntData &col_data)
 

Private Attributes

boost::shared_ptr< CommonDatacommonDataPtr
 

Detailed Description

Definition at line 67 of file RotatingFrameOperators.hpp.

Constructor & Destructor Documentation

◆ OpVolumeSideAssembleRowData()

OpRotatingFrameTools::OpVolumeSideAssembleRowData::OpVolumeSideAssembleRowData ( const std::string  field_name,
boost::shared_ptr< CommonData common_data_ptr 
)

Definition at line 255 of file RotatingFrameOperators.cpp.

258 commonDataPtr(common_data_ptr) {
259 sYmm = false;
260}
FormsIntegrators< DomainSideEleOp >::Assembly< USER_ASSEMBLE >::OpBase DomainSideEleOpAssembly
constexpr auto field_name
@ OPROWCOL
operator doWork is executed on FE rows &columns

Member Function Documentation

◆ doWork()

MoFEMErrorCode OpRotatingFrameTools::OpVolumeSideAssembleRowData::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
EntitiesFieldData::EntData row_data,
EntitiesFieldData::EntData col_data 
)

Definition at line 267 of file RotatingFrameOperators.cpp.

270 {
272 if (row_type != MBVERTEX)
274
275 const size_t nb_gauss_pts = getGaussPts().size2();
276 const size_t row_nb_dofs = commonDataPtr->faceRowData->getIndices().size();
277 const size_t col_nb_dofs = col_data.getIndices().size();
278
279 if (row_nb_dofs && col_nb_dofs) {
280
281 auto t_coords = getFTensor1CoordsAtGaussPts();
282 auto t_w = getFTensor0IntegrationWeight();
283 auto t_row_base = commonDataPtr->faceRowData->getFTensor0N();
284
285 auto t_row_diff_base = commonDataPtr->faceRowData->getFTensor1DiffN<3>();
286 size_t nb_face_functions = commonDataPtr->faceRowData->getN().size2();
287 auto t_normal = getFTensor1Normal();
288 const double ls = sqrt(t_normal(i) * t_normal(i));
289 t_normal(i) /= ls;
290
291 locMat.resize(row_nb_dofs, col_nb_dofs, false);
292 locMat.clear();
293
294 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
295 // we are actually integrating on a triangle here, therefore 0.5
296 // const double alpha_test = t_w * getMeasure() * (*cache).density;
297 const double alpha = t_w * commonDataPtr->mMeasure * (*cache).density;
298 // const double alpha = t_w * 0.5 * (*cache).density;
300 omg(i) = (*cache).Omega(i, j) * t_coords(j);
301
302 if (getNormalsAtGaussPts().size1() == nb_gauss_pts) {
303 VectorDouble n = getNormalsAtGaussPts(gg);
304 auto t_n = getFTensor1FromPtr<3>(&*n.data().begin());
305 t_normal(i) = t_n(i) / sqrt(t_n(j) * t_n(j));
306 }
307
308 size_t rr = 0;
309 for (; rr != row_nb_dofs / 3; ++rr) {
310 auto t_mat = getFTensor2FromArray<3, 3, 3>(locMat, 3 * rr);
311
312 const double beta = alpha * t_row_base;
313 auto t_col_diff_base = col_data.getFTensor1DiffN<3>(gg, 0);
314
315 for (size_t cc = 0; cc != col_nb_dofs / 3; ++cc) {
316
317 t_mat(i, k) += beta * kronecker_delta(i, k) *
318 (omg(j) * t_col_diff_base(j)) * (omg(l) * t_normal(l));
319
320 ++t_col_diff_base;
321 ++t_mat;
322 }
323
324 ++t_row_base;
325 }
326 for (; rr < nb_face_functions; ++rr)
327 ++t_row_base;
328
329 ++t_coords;
330 ++t_w;
331 }
332 // CHKERR iNtegrate(row_data, col_data);
333 CHKERR aSsemble(*commonDataPtr->faceRowData, col_data, false);
334 }
336}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
FTensor::Index< 'n', SPACE_DIM > n
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
Tensor2_Expr< Kronecker_Delta< T >, T, Dim0, Dim1, i, j > kronecker_delta(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
Rank 2.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
const VectorInt & getIndices() const
Get global indices of dofs on entity.

◆ iNtegrate()

MoFEMErrorCode OpRotatingFrameTools::OpVolumeSideAssembleRowData::iNtegrate ( EntData row_data,
EntData col_data 
)

Definition at line 262 of file RotatingFrameOperators.cpp.

263 {
264 return 0;
265}

Member Data Documentation

◆ commonDataPtr

boost::shared_ptr<CommonData> OpRotatingFrameTools::OpVolumeSideAssembleRowData::commonDataPtr
private

Definition at line 77 of file RotatingFrameOperators.hpp.


The documentation for this struct was generated from the following files: