v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
MoFEM::OpGetCoordsAndNormalsOnPrism Struct Reference

calculate normals at Gauss points of triangle element More...

#include <src/finite_elements/DataOperators.hpp>

Inheritance diagram for MoFEM::OpGetCoordsAndNormalsOnPrism:
[legend]
Collaboration diagram for MoFEM::OpGetCoordsAndNormalsOnPrism:
[legend]

Public Member Functions

 OpGetCoordsAndNormalsOnPrism (MatrixDouble &coords_at_gaussptf3, MatrixDouble &normals_at_gaussptf3, MatrixDouble &tangent1_at_gaussptf3, MatrixDouble &tangent2_at_gaussptf3, MatrixDouble &coords_at_gaussptf4, MatrixDouble &normals_at_gaussptf4, MatrixDouble &tangent1_at_gaussptf4, MatrixDouble &tangent2_at_gaussptf4)
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
MoFEMErrorCode calculateNormals ()
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Public Attributes

MatrixDoublecOords_at_GaussPtF3
 
MatrixDoublenOrmals_at_GaussPtF3
 
MatrixDoubletAngent1_at_GaussPtF3
 
MatrixDoubletAngent2_at_GaussPtF3
 
MatrixDoublecOords_at_GaussPtF4
 
MatrixDoublenOrmals_at_GaussPtF4
 
MatrixDoubletAngent1_at_GaussPtF4
 
MatrixDoubletAngent2_at_GaussPtF4
 
MatrixDouble sPin
 
- Public Attributes inherited from MoFEM::DataOperator
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 

Additional Inherited Members

- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 

Detailed Description

calculate normals at Gauss points of triangle element

Definition at line 378 of file DataOperators.hpp.

Constructor & Destructor Documentation

◆ OpGetCoordsAndNormalsOnPrism()

MoFEM::OpGetCoordsAndNormalsOnPrism::OpGetCoordsAndNormalsOnPrism ( MatrixDouble coords_at_gaussptf3,
MatrixDouble normals_at_gaussptf3,
MatrixDouble tangent1_at_gaussptf3,
MatrixDouble tangent2_at_gaussptf3,
MatrixDouble coords_at_gaussptf4,
MatrixDouble normals_at_gaussptf4,
MatrixDouble tangent1_at_gaussptf4,
MatrixDouble tangent2_at_gaussptf4 
)
inline

Definition at line 389 of file DataOperators.hpp.

394 : cOords_at_GaussPtF3(coords_at_gaussptf3),
395 nOrmals_at_GaussPtF3(normals_at_gaussptf3),
396 tAngent1_at_GaussPtF3(tangent1_at_gaussptf3),
397 tAngent2_at_GaussPtF3(tangent2_at_gaussptf3),
398 cOords_at_GaussPtF4(coords_at_gaussptf4),
399 nOrmals_at_GaussPtF4(normals_at_gaussptf4),
400 tAngent1_at_GaussPtF4(tangent1_at_gaussptf4),
401 tAngent2_at_GaussPtF4(tangent2_at_gaussptf4) {}

Member Function Documentation

◆ calculateNormals()

MoFEMErrorCode MoFEM::OpGetCoordsAndNormalsOnPrism::calculateNormals ( )

Definition at line 467 of file DataOperators.cpp.

467 {
469
470 sPin.resize(3, 3);
471 sPin.clear();
472 nOrmals_at_GaussPtF3.resize(tAngent1_at_GaussPtF3.size1(), 3, false);
473 for (unsigned int gg = 0; gg < tAngent1_at_GaussPtF3.size1(); ++gg) {
474 ierr = Spin(&*sPin.data().begin(), &tAngent1_at_GaussPtF3(gg, 0));
475 CHKERRG(ierr);
476 cblas_dgemv(CblasRowMajor, CblasNoTrans, 3, 3, 1., &*sPin.data().begin(), 3,
477 &tAngent2_at_GaussPtF3(gg, 0), 1, 0.,
478 &nOrmals_at_GaussPtF3(gg, 0), 1);
479 }
480 sPin.clear();
481 nOrmals_at_GaussPtF4.resize(tAngent1_at_GaussPtF4.size1(), 3, false);
482 for (unsigned int gg = 0; gg < tAngent1_at_GaussPtF4.size1(); ++gg) {
483 ierr = Spin(&*sPin.data().begin(), &tAngent1_at_GaussPtF4(gg, 0));
484 CHKERRG(ierr);
485 cblas_dgemv(CblasRowMajor, CblasNoTrans, 3, 3, 1., &*sPin.data().begin(), 3,
486 &tAngent2_at_GaussPtF4(gg, 0), 1, 0.,
487 &nOrmals_at_GaussPtF4(gg, 0), 1);
488 }
490}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
PetscErrorCode Spin(double *spinOmega, double *vecOmega)
calculate spin matrix from vector
Definition: fem_tools.c:546
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76

◆ doWork()

MoFEMErrorCode MoFEM::OpGetCoordsAndNormalsOnPrism::doWork ( int  side,
EntityType  type,
EntitiesFieldData::EntData data 
)
virtual

Operator for linear form, usually to calculate values on right hand side.

Reimplemented from MoFEM::DataOperator.

Definition at line 383 of file DataOperators.cpp.

384 {
386
387 if (data.getFieldData().size() == 0)
389 const int valid_edges3[] = {1, 1, 1, 0, 0, 0, 0, 0, 0};
390 const int valid_faces3[] = {0, 0, 0, 1, 0, 0, 0, 0, 0};
391 const int valid_edges4[] = {0, 0, 0, 0, 0, 0, 1, 1, 1};
392 const int valid_faces4[] = {0, 0, 0, 0, 1, 0, 0, 0, 0};
393
394 if (type == MBEDGE) {
395 if (!valid_edges3[side] || valid_edges4[side])
397 } else if (type == MBTRI) {
398 if (!valid_faces3[side] || valid_faces4[side])
400 }
401
402 switch (type) {
403 case MBVERTEX: {
404 for (unsigned int gg = 0; gg < data.getN().size1(); ++gg) {
405 for (int dd = 0; dd < 3; dd++) {
406 cOords_at_GaussPtF3(gg, dd) =
407 cblas_ddot(3, &data.getN(gg)[0], 1, &data.getFieldData()[dd], 3);
408 tAngent1_at_GaussPtF3(gg, dd) = cblas_ddot(
409 3, &data.getDiffN()(gg, 0), 2, &data.getFieldData()[dd], 3);
410 tAngent2_at_GaussPtF3(gg, dd) = cblas_ddot(
411 3, &data.getDiffN()(gg, 1), 2, &data.getFieldData()[dd], 3);
412 cOords_at_GaussPtF4(gg, dd) = cblas_ddot(
413 3, &data.getN(gg)[0], 1, &data.getFieldData()[9 + dd], 3);
414 tAngent1_at_GaussPtF4(gg, dd) = cblas_ddot(
415 3, &data.getDiffN()(gg, 6 + 0), 2, &data.getFieldData()[9 + dd], 3);
416 tAngent2_at_GaussPtF4(gg, dd) = cblas_ddot(
417 3, &data.getDiffN()(gg, 6 + 1), 2, &data.getFieldData()[9 + dd], 3);
418 }
419 }
420 } break;
421 case MBEDGE:
422 case MBTRI: {
423 if (2 * data.getN().size2() != data.getDiffN().size2()) {
424 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
425 }
426 unsigned int nb_dofs = data.getFieldData().size();
427 if (nb_dofs % 3 != 0) {
428 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
429 }
430 if (nb_dofs > 3 * data.getN().size2()) {
431 SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
432 "data inconsistency, side %d type %d", side, type);
433 }
434 for (unsigned int gg = 0; gg < data.getN().size1(); ++gg) {
435 for (int dd = 0; dd < 3; dd++) {
436 if ((type == MBTRI && valid_faces3[side]) ||
437 (type == MBEDGE && valid_edges3[side])) {
438 cOords_at_GaussPtF3(gg, dd) += cblas_ddot(
439 nb_dofs / 3, &data.getN(gg)[0], 1, &data.getFieldData()[dd], 3);
440 tAngent1_at_GaussPtF3(gg, dd) +=
441 cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 0), 2,
442 &data.getFieldData()[dd], 3);
443 tAngent2_at_GaussPtF3(gg, dd) +=
444 cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 1), 2,
445 &data.getFieldData()[dd], 3);
446 } else if ((type == MBTRI && valid_faces4[side]) ||
447 (type == MBEDGE && valid_edges4[side])) {
448 cOords_at_GaussPtF4(gg, dd) += cblas_ddot(
449 nb_dofs / 3, &data.getN(gg)[0], 1, &data.getFieldData()[dd], 3);
450 tAngent1_at_GaussPtF4(gg, dd) +=
451 cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 0), 2,
452 &data.getFieldData()[dd], 3);
453 tAngent2_at_GaussPtF4(gg, dd) +=
454 cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 1), 2,
455 &data.getFieldData()[dd], 3);
456 }
457 }
458 }
459 } break;
460 default:
461 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
462 }
463
465}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33

Member Data Documentation

◆ cOords_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::cOords_at_GaussPtF3

Definition at line 380 of file DataOperators.hpp.

◆ cOords_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::cOords_at_GaussPtF4

Definition at line 384 of file DataOperators.hpp.

◆ nOrmals_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::nOrmals_at_GaussPtF3

Definition at line 381 of file DataOperators.hpp.

◆ nOrmals_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::nOrmals_at_GaussPtF4

Definition at line 385 of file DataOperators.hpp.

◆ sPin

MatrixDouble MoFEM::OpGetCoordsAndNormalsOnPrism::sPin

Definition at line 403 of file DataOperators.hpp.

◆ tAngent1_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent1_at_GaussPtF3

Definition at line 382 of file DataOperators.hpp.

◆ tAngent1_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent1_at_GaussPtF4

Definition at line 386 of file DataOperators.hpp.

◆ tAngent2_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent2_at_GaussPtF3

Definition at line 383 of file DataOperators.hpp.

◆ tAngent2_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent2_at_GaussPtF4

Definition at line 387 of file DataOperators.hpp.


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