v0.13.2
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
BoneRemodeling::OpCalulatefRhoAtGaussPts Struct Reference

Assemble local vector containing density data. More...

#include <users_modules/bone_remodelling/src/DensityMaps.hpp>

Inheritance diagram for BoneRemodeling::OpCalulatefRhoAtGaussPts:
[legend]
Collaboration diagram for BoneRemodeling::OpCalulatefRhoAtGaussPts:
[legend]

Public Member Functions

 OpCalulatefRhoAtGaussPts (const string &row_field, CommonData &common_data, DataFromMetaIO &data_from_meta_io, SurfaceKDTree &surface_distance)
 
PetscErrorCode doWork (int row_side, EntityType row_type, DataForcesAndSourcesCore::EntData &row_data)
 
- 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)
 User call this function to loop over elements on the side of face. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is 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 call this function to loop over parent elements. This function calls finite element with is 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 call this function to loop over parent elements. This function calls finite element with is 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

CommonDatacommonData
 
DataFromMetaIOdataFromMetaIo
 
SurfaceKDTreesurfaceDistance
 
- 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...
 
- 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

Assemble local vector containing density data.

Definition at line 769 of file DensityMaps.hpp.

Constructor & Destructor Documentation

◆ OpCalulatefRhoAtGaussPts()

BoneRemodeling::OpCalulatefRhoAtGaussPts::OpCalulatefRhoAtGaussPts ( const string &  row_field,
CommonData common_data,
DataFromMetaIO data_from_meta_io,
SurfaceKDTree surface_distance 
)
inline

Definition at line 775 of file DensityMaps.hpp.

780 :
782 row_field,row_field,UserDataOperator::OPROW
783 ),
784 commonData(common_data),
785 dataFromMetaIo(data_from_meta_io),
786 surfaceDistance(surface_distance) {
787 }

Member Function Documentation

◆ doWork()

PetscErrorCode BoneRemodeling::OpCalulatefRhoAtGaussPts::doWork ( int  row_side,
EntityType  row_type,
DataForcesAndSourcesCore::EntData &  row_data 
)
inline

< reference vectors of coordinates

< reference vector of indices

< reference vector of data (colors)

< vector of the densities

< vector of reference distances of the points

Definition at line 790 of file DensityMaps.hpp.

792 {
793 PetscFunctionBegin;
794
795 if(row_type!=MBVERTEX) PetscFunctionReturn(0);
796
797 int nb_rows = row_data.getIndices().size();
798 if(nb_rows == 0) PetscFunctionReturn(0);
799
800 ////////////////////////////////////////////////////////////////////////////////////////////
801 vector<int> vecIx; ///< reference vectors of coordinates
802 vector<int> vecIy;
803 vector<int> vecIz;
804 vector<int> vecIdx; ///< reference vector of indices
805
806 vector<double> vecData; ///< reference vector of data (colors)
807 vector<double> vecDensity; ///< vector of the densities
808 vector<double> vecDist; ///< vector of reference distances of the points
809 ///////////////////////////////////////////////////////////////////////////////////////////////////
810
811 auto t0_dist = getFTensor0FromVec(*commonData.distanceValue);
812 auto t1_grad = getFTensor1FromMat<3>(*commonData.distanceGradient);
814
816 &getCoordsAtGaussPts()(0,0),
817 &getCoordsAtGaussPts()(0,1),
818 &getCoordsAtGaussPts()(0,2),
819 3
820 );
821
823 dataFromMetaIo.metaFile.ElementSize()[0] * dataFromMetaIo.sCale,
824 dataFromMetaIo.metaFile.ElementSize()[1] * dataFromMetaIo.sCale,
826 );
827 t1_dx(i) *= dataFromMetaIo.cUbe_size;
828 double cube_size = sqrt(t1_dx(i)*t1_dx(i))*dataFromMetaIo.dIst_tol;
829
830 const int nb_gauss_pts = row_data.getN().size1();
831 commonData.rhoAtGaussPts.resize(nb_gauss_pts,false);
832 for(int gg = 0;gg<nb_gauss_pts;gg++) {
833
834 if(t0_dist<cube_size) {
835 const double nrm2_grad = sqrt(t1_grad(i)*t1_grad(i));
837 t1_delta(i) = t1_grad(i)/nrm2_grad;
838 t1_xyz(i) -= (t0_dist-cube_size)*t1_delta(i);
839 }
840
841
842 // //double cube_size = sqrt(dx*dx+dy*dy+dz*dz)*dataFromMetaIo.dIst_tol; //FIXME
843 // double xc_dm = (2 * dataFromMetaIo.cUbe_size - 1) * dataFromMetaIo.metaFile.ElementSize()[0];
844 // double yc_dm = (2 * dataFromMetaIo.cUbe_size - 1) * dataFromMetaIo.metaFile.ElementSize()[1];
845 // double zc_dm = (2 * dataFromMetaIo.cUbe_size - 1) * dataFromMetaIo.metaFile.ElementSize()[2];
846 // double cube_size = 0.5 * sqrt(xc_dm * xc_dm + yc_dm * yc_dm + zc_dm * zc_dm) * dataFromMetaIo.dIst_tol; // real cube size /2
847 //
848 // if(dist<cube_size) {
849 // // Shift cube from boundary
850 // x += (dist-cube_size)*delta_x/dist;
851 // y += (dist-cube_size)*delta_y/dist;
852 // z += (dist-cube_size)*delta_z/dist;
853 // }
854
855 //Calculate density like in OpMassCalculation%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
857 t1_xyz(0),t1_xyz(1),t1_xyz(2),
858 t1_dx(0),t1_dx(1),t1_dx(2),
859 vecIx,vecIy,vecIz,vecIdx, vecDist
860 ); CHKERRQ(ierr);
861
862 ierr = dataFromMetaIo.getDataForGiveIndex(vecIdx,vecData); CHKERRQ(ierr);
863 ierr = dataFromMetaIo.mapDensityLinear(vecData,vecDensity); CHKERRQ(ierr);
864
865 double data = 0.0;
866 switch (dataFromMetaIo.fIlter) {
867 case 0:
868 data = dataFromMetaIo.mapGaussSmooth(vecDensity,vecDist);
869 break;
870 case 1:
871 data = dataFromMetaIo.getAverage(vecDensity);
872 break;
873 case 2:
874 data = dataFromMetaIo.medianFilter(vecDensity);
875
876 default:
877 data = dataFromMetaIo.getAverage(vecDensity);
878 }
879
880 commonData.rhoAtGaussPts[gg] = data;//pow(x,2)*y; //data; // not multiply this by vol
881
882 ++t1_xyz;
883 ++t0_dist;
884 ++t1_grad;
885
886 }
887
888 PetscFunctionReturn(0);
889 }
static PetscErrorCode ierr
FTensor::Index< 'i', SPACE_DIM > i
VectorDouble rhoAtGaussPts
Values of density at integration Pts.
boost::shared_ptr< MatrixDouble > distanceGradient
boost::shared_ptr< VectorDouble > distanceValue
PetscErrorCode getDataForGiveIndex(const vector< int > &vec_idx, vector< double > &vec_data)
double medianFilter(const vector< double > &vec_density)
function returns median of data from given vector
PetscErrorCode mapDensityLinear(double *paramDensity, const vector< double > &vec_data, vector< double > &vec_density)
Transform data to densities.
PetscErrorCode getInidcesForGivenCoordinate(const double x, const double y, const double z, const double dx, const double dy, const double dz, vector< int > &vec_ix, vector< int > &vec_iy, vector< int > &vec_iz, vector< int > &vec_idx, vector< double > &vec_dist)
double getAverage(const vector< double > &vec_data)
Averages data within given vector.
double mapGaussSmooth(const double sigma, const vector< double > &vec_density, const vector< double > &vec_dist)
function smooths values in a given vector depending on distance from center (vec_dist)
MatrixDouble & getCoordsAtGaussPts()
Gauss points and weight, matrix (nb. of points x 3)

Member Data Documentation

◆ commonData

CommonData& BoneRemodeling::OpCalulatefRhoAtGaussPts::commonData

Definition at line 771 of file DensityMaps.hpp.

◆ dataFromMetaIo

DataFromMetaIO& BoneRemodeling::OpCalulatefRhoAtGaussPts::dataFromMetaIo

Definition at line 772 of file DensityMaps.hpp.

◆ surfaceDistance

SurfaceKDTree& BoneRemodeling::OpCalulatefRhoAtGaussPts::surfaceDistance

Definition at line 773 of file DensityMaps.hpp.


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