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

#include <users_modules/homogenisation/src/NitschePeriodicMethod.hpp>

Inheritance diagram for PeriodicNitscheConstrains::OpLhsPeriodicNormal:
[legend]
Collaboration diagram for PeriodicNitscheConstrains::OpLhsPeriodicNormal:
[legend]

Public Member Functions

 OpLhsPeriodicNormal (const string field_name, NitscheMethod::BlockData &nitsche_block_data, NitscheMethod::CommonData &nitsche_common_data, NonlinearElasticElement::BlockData &data, NonlinearElasticElement::CommonData &common_data, CommonData &periodic_common_data, bool field_disp=true)
 
PetscErrorCode doWork (int row_side, EntityType row_type, DataForcesAndSourcesCore::EntData &row_data)
 
- Public Member Functions inherited from NitscheMethod::OpCommon
 OpCommon (const std::string field_name, BlockData &nitsche_block_data, CommonData &nitsche_common_data, NonlinearElasticElement::BlockData &data, NonlinearElasticElement::CommonData &common_data, bool field_disp, const char type)
 
virtual MoFEMErrorCode calculateP (int gg, int fgg, int ff)
 
MoFEMErrorCode getJac (DataForcesAndSourcesCore::EntData &data, int gg, MatrixDouble &jac)
 
MoFEMErrorCode getTractionVariance (int gg, int fgg, int ff, MatrixDouble &jac, MatrixDouble &trac)
 
- Public Member Functions inherited from NitscheMethod::OpBasicCommon
 OpBasicCommon (const std::string field_name, BlockData &nitsche_block_data, CommonData &nitsche_common_data, bool field_disp, const char type)
 
virtual MoFEMErrorCode getFaceRadius (int ff)
 
virtual MoFEMErrorCode getGammaH (double gamma, int gg)
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
doublegetVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
VolumeElementForcesAndSourcesCoregetVolumeFE () const
 return pointer to Generic Volume Finite Element object More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
std::string getFEName () const
 Get name of the element. More...
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. More...
 
- 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

CommonDataperiodicCommonData
 
MatrixDouble kMatrix0
 
MatrixDouble kMatrix1
 
MatrixDouble kMatrixOff
 
- Public Attributes inherited from NitscheMethod::OpCommon
NonlinearElasticElement::BlockDatadAta
 
NonlinearElasticElement::CommonDatacommonData
 
VectorDouble dIsp
 
VectorDouble tRaction
 
MatrixDouble jAc_row
 
MatrixDouble jAc_col
 
MatrixDouble tRac_v
 
MatrixDouble tRac_u
 
- Public Attributes inherited from NitscheMethod::OpBasicCommon
BlockDatanitscheBlockData
 
CommonDatanitscheCommonData
 
bool fieldDisp
 
double faceRadius
 
double gammaH
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- 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::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 1 << 3 ,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
using AdjCache = std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > > >
 
- 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)>
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Definition at line 740 of file NitschePeriodicMethod.hpp.

Constructor & Destructor Documentation

◆ OpLhsPeriodicNormal()

PeriodicNitscheConstrains::OpLhsPeriodicNormal::OpLhsPeriodicNormal ( const string  field_name,
NitscheMethod::BlockData nitsche_block_data,
NitscheMethod::CommonData nitsche_common_data,
NonlinearElasticElement::BlockData data,
NonlinearElasticElement::CommonData common_data,
CommonData periodic_common_data,
bool  field_disp = true 
)
inline

Definition at line 744 of file NitschePeriodicMethod.hpp.

752 :
753 OpCommon(
755 nitsche_block_data,
756 nitsche_common_data,
757 data,
758 common_data,
759 field_disp,
760 UserDataOperator::OPROW
761 ),
762 periodicCommonData(periodic_common_data) {
763 }
constexpr auto field_name
OpCommon(const std::string field_name, BlockData &nitsche_block_data, CommonData &nitsche_common_data, NonlinearElasticElement::BlockData &data, NonlinearElasticElement::CommonData &common_data, bool field_disp, const char type)

Member Function Documentation

◆ doWork()

PetscErrorCode PeriodicNitscheConstrains::OpLhsPeriodicNormal::doWork ( int  row_side,
EntityType  row_type,
DataForcesAndSourcesCore::EntData &  row_data 
)
inline

Definition at line 768 of file NitschePeriodicMethod.hpp.

770 {
771 PetscFunctionBegin;
772
773
774 if(dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) == dAta.tEts.end()) {
775 PetscFunctionReturn(0);
776 }
777 if(row_data.getIndices().size()==0) PetscFunctionReturn(0);
778 int nb_dofs_row = row_data.getIndices().size();
779 const double gamma = nitscheBlockData.gamma;
780 const double phi = nitscheBlockData.phi;
781 double eps = periodicCommonData.eps;
782
783 try {
784
785 int gg = 0;
786 for(int ff = 0;ff<4;ff++) {
787 if(nitscheCommonData.facesFePtr[ff]==NULL) continue;
788 int nb_face_gauss_pts = nitscheCommonData.faceGaussPts[ff].size2();
789 vector<vector<MatrixDouble > > &P = nitscheCommonData.P;
790 P.resize(4);
791 P[ff].resize(nb_face_gauss_pts);
792 for(int fgg = 0;fgg<nb_face_gauss_pts;fgg++) {
793 P[ff][fgg].resize(3,3);
794 P[ff][fgg].clear();
795 for(int dd = 0;dd<3;dd++) {
796 P[ff][fgg](dd,dd) = 1;
797 }
798 }
799 ierr = getFaceRadius(ff); CHKERRQ(ierr);
800 for(int fgg = 0;fgg<nb_face_gauss_pts;fgg++,gg++) {
801 ierr = getGammaH(gamma,gg); CHKERRQ(ierr);
802 double val = getGaussPts()(3,gg);
803 ierr = getJac(row_data,gg,jAc_row); CHKERRQ(ierr);
804 ierr = getTractionVariance(gg,fgg,ff,jAc_row,tRac_v); CHKERRQ(ierr);
806 3,ublas::shallow_array_adaptor<double>(3,&nitscheCommonData.faceNormals[ff](fgg,0))
807 );
808 double area = cblas_dnrm2(3,&normal[0],1);
809 const double alpha = 0.5; // mortar parameter
810
811 CommonData::Container::nth_index<0>::type::iterator it,hi_it;
812 it = periodicCommonData.facesContainer.get<0>().lower_bound(gg);
813 hi_it = periodicCommonData.facesContainer.get<0>().upper_bound(gg);
814 for(;it!=hi_it;it++) {
815
816 const ublas::vector<int> &indices = it->iNdices;
817 const VectorDouble &shape = it->shapeFunctions;
818
819 int nb_dofs_col = indices.size();
820 if(indices.size()==0) {
821 continue;
822 }
823 kMatrix0.resize(nb_dofs_row,nb_dofs_col,false);
824 kMatrix0.clear();
825 for(int dd1 = 0;dd1<nb_dofs_row/3;dd1++) {
826 double n_row = row_data.getN()(gg,dd1);
827 for(int dd2 = 0;dd2<nb_dofs_col/3;dd2++) {
828 double n_col = shape[dd2];
829 for(int dd3 = 0;dd3<3;dd3++) {
830 kMatrix0(3*dd1+dd3,3*dd2+dd3) -= (1-eps)*val*n_row*n_col*area;
831 }
832 }
833 }
834 kMatrix1.resize(nb_dofs_row,nb_dofs_col,false);
835 kMatrix1.clear();
836 for(int dd1 = 0;dd1<nb_dofs_row;dd1++) {
837 for(int dd2 = 0;dd2<nb_dofs_col/3;dd2++) {
838 double n_col = shape[dd2];
839 for(int dd3 = 0;dd3<3;dd3++) {
840 kMatrix1(dd1,3*dd2+dd3) += (1-eps)*phi*val*tRac_v(dd3,dd1)*n_col;
841 }
842 }
843 }
844
845 kMatrix0 /= gammaH;
847 getFEMethod()->snes_B,
848 nb_dofs_row,
849 &row_data.getIndices()[0],
850 nb_dofs_col,
851 &indices[0],
852 &kMatrix0(0,0),
853 ADD_VALUES
854 ); CHKERRQ(ierr);
855
857 getFEMethod()->snes_B,
858 nb_dofs_row,
859 &row_data.getIndices()[0],
860 nb_dofs_col,
861 &indices[0],
862 &kMatrix1(0,0),
863 ADD_VALUES
864 ); CHKERRQ(ierr);
865
866 kMatrix1.resize(nb_dofs_col,nb_dofs_row);
867 kMatrix1.clear();
868 for(int dd1 = 0;dd1<nb_dofs_col/3;dd1++) {
869 double n_col = shape[dd1];
870 for(int dd2 = 0;dd2<3;dd2++) {
871 for(int dd3 = 0;dd3<nb_dofs_row;dd3++) {
872 kMatrix1(3*dd1+dd2,dd3) += (1-eps)*val*n_col*tRac_v(dd2,dd3);
873 }
874 }
875 }
876
877 kMatrix1 *= alpha;
879 getFEMethod()->snes_B,
880 nb_dofs_col,
881 &indices[0],
882 nb_dofs_row,
883 &row_data.getIndices()[0],
884 &kMatrix1(0,0),
885 ADD_VALUES
886 ); CHKERRQ(ierr);
887
888 }
889
890 MatrixDouble jac;
891 MatrixDouble trac;
892
893 it = periodicCommonData.volumesContainer.get<0>().lower_bound(gg);
894 hi_it = periodicCommonData.volumesContainer.get<0>().upper_bound(gg);
895 for(;it!=hi_it;it++) {
896
897 int nb = it->iNdices.size();
898 jac.resize(9,nb,false);
899 jac.clear();
900 const MatrixDouble diff_n = it->diffShapeFunctions;
901 const MatrixDouble jac_stress = periodicCommonData.stressJacobian[gg];
902 for(int dd = 0;dd<nb/3;dd++) {
903 for(int rr = 0;rr<3;rr++) {
904 for(int ii = 0;ii<9;ii++) {
905 for(int jj = 0;jj<3;jj++) {
906 jac(ii,3*dd+rr) += jac_stress(ii,3*rr+jj)*diff_n(dd,jj);
907 }
908 }
909 }
910 }
911 trac.resize(3,jac.size2());
912 trac.clear();
913 for(unsigned int dd2 = 0;dd2<jac.size2();dd2++) {
914 for(int nn = 0;nn<3;nn++) {
915 trac(nn,dd2) = -cblas_ddot(3,&jac(3*nn,dd2),jac.size2(),&normal[0],1);
916 }
917 }
918
919 kMatrixOff.resize(nb_dofs_row,nb);
920 kMatrixOff.clear();
921 for(int dd1 = 0;dd1<nb_dofs_row/3;dd1++) {
922 double n_row = row_data.getN()(gg,dd1);
923 for(int dd2 = 0;dd2<3;dd2++) {
924 for(int dd3 = 0;dd3<nb;dd3++) {
925 kMatrixOff(3*dd1+dd2,dd3) += (1-eps)*val*n_row*trac(dd2,dd3);
926 }
927 }
928 }
929
930 kMatrixOff *= (1-alpha);
932 getFEMethod()->snes_B,
933 nb_dofs_row,
934 &row_data.getIndices()[0],
935 nb,
936 &it->iNdices[0],
937 &kMatrixOff(0,0),
938 ADD_VALUES
939 ); CHKERRQ(ierr);
940
941 }
942
943 }
944 }
945
946 } catch (const std::exception& ex) {
947 ostringstream ss;
948 ss << "throw in method: " << ex.what() << endl;
949 SETERRQ(PETSC_COMM_SELF,1,ss.str().c_str());
950 }
951
952 PetscFunctionReturn(0);
953 }
static PetscErrorCode ierr
static const double eps
static double phi
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
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
FTensor::Tensor2< double *, 3, 3 > & getJac()
get element Jacobian
double gamma
Penalty term, see .
double phi
Nitsche method parameter, see .
std::vector< boost::shared_ptr< const NumeredEntFiniteElement > > facesFePtr
std::vector< MatrixDouble > faceNormals
std::vector< std::vector< MatrixDouble > > P
projection matrix
std::vector< MatrixDouble > faceGaussPts
virtual MoFEMErrorCode getFaceRadius(int ff)
virtual MoFEMErrorCode getGammaH(double gamma, int gg)
NonlinearElasticElement::BlockData & dAta
MoFEMErrorCode getTractionVariance(int gg, int fgg, int ff, MatrixDouble &jac, MatrixDouble &trac)
Range tEts
constrains elements in block set

Member Data Documentation

◆ kMatrix0

MatrixDouble PeriodicNitscheConstrains::OpLhsPeriodicNormal::kMatrix0

Definition at line 765 of file NitschePeriodicMethod.hpp.

◆ kMatrix1

MatrixDouble PeriodicNitscheConstrains::OpLhsPeriodicNormal::kMatrix1

Definition at line 765 of file NitschePeriodicMethod.hpp.

◆ kMatrixOff

MatrixDouble PeriodicNitscheConstrains::OpLhsPeriodicNormal::kMatrixOff

Definition at line 766 of file NitschePeriodicMethod.hpp.

◆ periodicCommonData

CommonData& PeriodicNitscheConstrains::OpLhsPeriodicNormal::periodicCommonData

Definition at line 742 of file NitschePeriodicMethod.hpp.


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