v0.15.5
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
MoFEM::Projection10NodeCoordsFromFineSurfaceOnField Struct Reference

Project edge mid-nodes onto a fine surface mesh using closest points. More...

#include "src/approximation/Projection10NodeCoordsOnField.hpp"

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

Public Member Functions

 Projection10NodeCoordsFromFineSurfaceOnField (Interface &m_field, moab::Interface &fine_moab, const Range &fine_tris, std::string field_name, double planar_angle_deg=5.0, double max_disp_factor=0.5, bool debug_surfaces=false, int verb=0)
 
MoFEMErrorCode preProcess ()
 Pre-processing function executed at loop initialization.
 
MoFEMErrorCode operator() ()
 Main operator function executed for each loop iteration.
 
MoFEMErrorCode postProcess ()
 Post-processing function executed at loop completion.
 
- Public Member Functions inherited from MoFEM::EntityMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 EntityMethod ()=default
 Default constructor.
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 Default constructor.
 
virtual ~BasicMethod ()=default
 Virtual destructor.
 
int getNinTheLoop () const
 Get current loop iteration index.
 
int getLoopSize () const
 Get total loop size.
 
auto getLoHiFERank () const
 Get processor rank range for finite element iteration.
 
auto getLoFERank () const
 Get lower processor rank for finite element iteration.
 
auto getHiFERank () const
 Get upper processor rank for finite element iteration.
 
unsigned int getFieldBitNumber (std::string field_name) const
 Get bit number for a specific field by name.
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from another BasicMethod instance.
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak pointer object.
 
- Public Member Functions inherited from MoFEM::KspMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 KspMethod ()
 Default constructor.
 
virtual ~KspMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 Copy data from another KSP method.
 
- Public Member Functions inherited from MoFEM::PetscData
 PetscData ()
 Default constructor.
 
virtual ~PetscData ()=default
 Virtual destructor.
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 Copy PETSc data from another instance.
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference to pointer of interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface.
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface.
 
virtual ~UnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::SnesMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 SnesMethod ()
 Default constructor.
 
virtual ~SnesMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy SNES data from another instance.
 
- Public Member Functions inherited from MoFEM::TSMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 TSMethod ()
 Default constructor.
 
virtual ~TSMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data from another instance.
 
- Public Member Functions inherited from MoFEM::TaoMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 TaoMethod ()
 Default constructor.
 
virtual ~TaoMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyTao (const TaoMethod &tao)
 Copy TAO data from another instance.
 

Protected Member Functions

 FTENSOR_INDEX (3, i)
 
 FTENSOR_INDEX (3, j)
 

Protected Attributes

InterfacemField
 
moab::Interface & fineMoab
 
Range fineTris
 
std::string fieldName
 
int vErbose
 
AdaptiveKDTree fineTree
 
EntityHandle fineRoot
 
AdaptiveKDTree fineFeatureTree
 
EntityHandle fineFeatureRoot
 
Range fineFeatureTris
 
double planarAngleCos
 
double maxDispFactor
 
bool debugSurfaces
 
Range coarseSkinTris
 
Range coarseSkinEdges
 
VectorDouble coords
 
VectorDouble3 aveMidCoord
 
VectorDouble3 midNodeCoord
 
VectorDouble3 diffNodeCoord
 
VectorDouble3 dOf
 

Additional Inherited Members

- Public Types inherited from MoFEM::KspMethod
enum  KSPContext { CTX_SETFUNCTION , CTX_OPERATORS , CTX_KSPNONE }
 Context enumeration for KSP solver phases. More...
 
- Public Types inherited from MoFEM::PetscData
enum  DataContext {
  CTX_SET_NONE = 0 , CTX_SET_F = 1 << 0 , CTX_SET_A = 1 << 1 , CTX_SET_B = 1 << 2 ,
  CTX_SET_X = 1 << 3 , CTX_SET_DX = 1 << 4 , CTX_SET_X_T = 1 << 5 , CTX_SET_X_TT = 1 << 6 ,
  CTX_SET_TIME = 1 << 7
}
 Enumeration for data context flags. More...
 
using Switches = std::bitset< 8 >
 Bitset type for context switches.
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION , CTX_SNESSETJACOBIAN , CTX_SNESNONE }
 Context enumeration for SNES solver phases. More...
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION , CTX_TSSETRHSJACOBIAN , CTX_TSSETIFUNCTION , CTX_TSSETIJACOBIAN ,
  CTX_TSTSMONITORSET , CTX_TSNONE
}
 Context enumeration for TS solver phases. More...
 
- Public Types inherited from MoFEM::TaoMethod
enum  TAOContext { CTX_TAO_OBJECTIVE , CTX_TAO_GRADIENT , CTX_TAO_HESSIAN , CTX_TAO_NONE }
 Context enumeration for TAO solver phases. More...
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version.
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version.
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version.
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version.
 
- Public Attributes inherited from MoFEM::EntityMethod
boost::shared_ptr< FieldfieldPtr
 Shared pointer to field information.
 
boost::shared_ptr< FieldEntityentPtr
 Shared pointer to field entity data.
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 Current index of processed method in the loop.
 
int loopSize
 Total number of methods to process in the loop.
 
std::pair< int, int > loHiFERank
 Processor rank range for distributed finite element iteration.
 
int rAnk
 Current processor rank in MPI communicator.
 
int sIze
 Total number of processors in MPI communicator.
 
const RefEntity_multiIndexrefinedEntitiesPtr
 Pointer to container of refined MoFEM DOF entities.
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 Pointer to container of refined finite element entities.
 
const ProblemproblemPtr
 Raw pointer to current MoFEM problem instance.
 
const Field_multiIndexfieldsPtr
 Raw pointer to fields multi-index container.
 
const FieldEntity_multiIndexentitiesPtr
 Raw pointer to container of field entities.
 
const DofEntity_multiIndexdofsPtr
 Raw pointer to container of degree of freedom entities.
 
const FiniteElement_multiIndexfiniteElementsPtr
 Raw pointer to container of finite elements.
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 Raw pointer to container of finite element entities.
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 Raw pointer to container of adjacencies between DOFs and finite elements.
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing operations.
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing operations.
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for main operator execution.
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 Switch for vector assembly operations.
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 Switch for matrix assembly operations.
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 Weak pointer to cached entity data.
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Current KSP computation context.
 
KSP ksp
 PETSc KSP linear solver object.
 
Vec & ksp_f
 Reference to residual vector in KSP context.
 
Mat & ksp_A
 Reference to system matrix in KSP context.
 
Mat & ksp_B
 Reference to preconditioner matrix in KSP context.
 
- Public Attributes inherited from MoFEM::PetscData
Switches data_ctx
 Current data context switches.
 
Vec f
 PETSc residual vector.
 
Mat A
 PETSc Jacobian matrix.
 
Mat B
 PETSc preconditioner matrix.
 
Vec x
 PETSc solution vector.
 
Vec dx
 PETSc solution increment vector.
 
Vec x_t
 PETSc first time derivative vector.
 
Vec x_tt
 PETSc second time derivative vector.
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 Current SNES computation context.
 
SNES snes
 PETSc SNES nonlinear solver object.
 
Vec & snes_x
 Reference to current solution state vector.
 
Vec & snes_dx
 Reference to solution update/increment vector.
 
Vec & snes_f
 Reference to residual vector.
 
Mat & snes_A
 Reference to Jacobian matrix.
 
Mat & snes_B
 Reference to preconditioner of Jacobian matrix.
 
- Public Attributes inherited from MoFEM::TSMethod
TS ts
 PETSc time stepping solver object.
 
TSContext ts_ctx
 Current TS computation context.
 
PetscInt ts_step
 Current time step number.
 
PetscReal ts_a
 Shift parameter for U_t (see PETSc Time Solver documentation)
 
PetscReal ts_aa
 Shift parameter for U_tt (second time derivative)
 
PetscReal ts_t
 Current time value.
 
PetscReal ts_dt
 Current time step size.
 
Vec & ts_u
 Reference to current state vector U(t)
 
Vec & ts_u_t
 Reference to first time derivative of state vector dU/dt.
 
Vec & ts_u_tt
 Reference to second time derivative of state vector d²U/dt²
 
Vec & ts_F
 Reference to residual vector F(t,U,U_t,U_tt)
 
Mat & ts_A
 Reference to Jacobian matrix: dF/dU + a*dF/dU_t + aa*dF/dU_tt.
 
Mat & ts_B
 Reference to preconditioner matrix for ts_A.
 
- Public Attributes inherited from MoFEM::TaoMethod
TAOContext tao_ctx
 Current TAO computation context.
 
Tao tao
 PETSc TAO optimization solver object.
 
Vec & tao_x
 Reference to optimization variables vector.
 
Vec & tao_f
 Reference to gradient vector.
 
Mat & tao_A
 Reference to Hessian matrix.
 
Mat & tao_B
 Reference to preconditioner matrix for Hessian.
 
- Static Public Attributes inherited from MoFEM::PetscData
static constexpr Switches CtxSetNone = PetscData::Switches(CTX_SET_NONE)
 No data switch.
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 Residual vector switch.
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 Jacobian matrix switch.
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 Preconditioner matrix switch.
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 Solution vector switch.
 
static constexpr Switches CtxSetDX = PetscData::Switches(CTX_SET_DX)
 Solution increment switch.
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 First time derivative switch.
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 Second time derivative switch.
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 Time value switch.
 

Detailed Description

Project edge mid-nodes onto a fine surface mesh using closest points.

Definition at line 51 of file Projection10NodeCoordsOnField.hpp.

Constructor & Destructor Documentation

◆ Projection10NodeCoordsFromFineSurfaceOnField()

MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::Projection10NodeCoordsFromFineSurfaceOnField ( Interface m_field,
moab::Interface &  fine_moab,
const Range fine_tris,
std::string  field_name,
double  planar_angle_deg = 5.0,
double  max_disp_factor = 0.5,
bool  debug_surfaces = false,
int  verb = 0 
)

Definition at line 194 of file Projection10NodeCoordsOnField.cpp.

203 : mField(m_field), fineMoab(fine_moab), fineTris(fine_tris),
204 fieldName(field_name), vErbose(verb), fineTree(&fine_moab), fineRoot(0),
205 fineFeatureTree(&fine_moab), fineFeatureRoot(0),
206 planarAngleCos(std::cos(planar_angle_deg * 0.017453292519943295)),
207 maxDispFactor(max_disp_factor), debugSurfaces(debug_surfaces) {}
constexpr auto field_name

Member Function Documentation

◆ FTENSOR_INDEX() [1/2]

MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::FTENSOR_INDEX ( ,
i   
)
protected

◆ FTENSOR_INDEX() [2/2]

MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::FTENSOR_INDEX ( ,
j   
)
protected

◆ operator()()

MoFEMErrorCode MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::operator() ( )
virtual

Main operator function executed for each loop iteration.

This virtual function is called for every item (finite element, entity, etc.) in the loop. It is the core computation function typically used for:

  • Calculating element local matrices and vectors
  • Assembling contributions to global system
  • Element-level post-processing operations
Returns
Error code indicating success or failure

Reimplemented from MoFEM::BasicMethod.

Definition at line 302 of file Projection10NodeCoordsOnField.cpp.

302 {
304
305 if (entPtr == NULL) {
306 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
307 }
308 if (entPtr->getName() != fieldName)
310 if (entPtr->getEntType() != MBEDGE)
312
313 const EntityHandle edge = entPtr->getEnt();
314 if (coarseSkinEdges.find(edge) == coarseSkinEdges.end()) {
316 }
317
318 Range adj_tris;
319 CHKERR mField.get_moab().get_adjacencies(&edge, 1, 2, false, adj_tris,
320 moab::Interface::UNION);
321 adj_tris = intersect(adj_tris, coarseSkinTris);
322 if (adj_tris.size() != 2) {
323 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
324 "skin edge does not have exactly two adjacent triangles");
325 }
326
327 int num_nodes = 0;
328 const EntityHandle *conn = nullptr;
329 CHKERR mField.get_moab().get_connectivity(edge, conn, num_nodes, false);
330 if (num_nodes < 2) {
331 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
332 "edge connectivity has fewer than 2 nodes");
333 }
334
335 const auto base = entPtr->getApproxBase();
336 double edge_shape_function_val = 0.0;
337 switch (base) {
339 edge_shape_function_val = 0.25;
340 break;
342 edge_shape_function_val = 0.25 * LOBATTO_PHI0(0);
343 break;
345 double L[3];
346 CHKERR Legendre_polynomials(2, 0, NULL, L, NULL, 1);
347 edge_shape_function_val = 0.125 * LOBATTO_PHI0(0);
348 } break;
349 default:
350 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not yet implemented");
351 }
352
353 auto ent_data = entPtr->getEntFieldData();
354 auto t_ent_data = getFTensor1FromPtr<3>(&ent_data(0));
355 if (ent_data.size() != 3) {
356 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
357 "edge data size mismatch");
358 }
359 coords.resize(num_nodes * 3);
360 CHKERR mField.get_moab().get_coords(conn, num_nodes, &*coords.data().begin());
363 &coords[0], &coords[1], &coords[2]);
365 &coords[3], &coords[4], &coords[5]);
366 edge_vec(i) = t_n1(i) - t_n0(i);
367 const double h = edge_vec.l2();
368 if (h <= 0.0) {
369 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "zero-length edge");
370 }
371 aveMidCoord.resize(3);
372 midNodeCoord.resize(3);
373 auto t_ave_mid_coord = getFTensor1FromPtr<3>(&aveMidCoord(0));
374 auto t_mid_node_coord = getFTensor1FromPtr<3>(&midNodeCoord(0));
375 t_ave_mid_coord(i) = 0.5 * (t_n0(i) + t_n1(i));
376 t_mid_node_coord(i) =
377 t_ave_mid_coord(i) + edge_shape_function_val * t_ent_data(i);
378
379 auto tri_normal = [&](const EntityHandle tri,
380 FTensor::Tensor1<double, 3> &t_normal) {
382 CHKERR mField.getInterface<Tools>()->getTriNormal(tri, &t_normal(0));
383 t_normal.normalize();
385 };
386
389 CHKERR tri_normal(adj_tris[0], n0);
390 CHKERR tri_normal(adj_tris[1], n1);
391 const double dot = n0(i) * n1(i);
392
393 if (std::abs(dot) >= planarAngleCos) {
395 }
396
397 double closest_point[3] = {0.0, 0.0, 0.0};
398 EntityHandle closest_tri = 0;
399 const double query[3] = {midNodeCoord[0], midNodeCoord[1], midNodeCoord[2]};
400 CHKERR fineFeatureTree.closest_triangle(fineFeatureRoot, query, closest_point,
401 closest_tri);
402 if (!closest_tri) {
403 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
404 "no closest triangle found on fine mesh");
405 }
406
407 auto get_disp = [&](double _closest_point[3]) {
409 for (int ii = 0; ii < 3; ii++)
410 delta(ii) = _closest_point[ii] - midNodeCoord[ii];
411 return delta.l2();
412 };
413
414 double disp = get_disp(closest_point);
415
416 if (disp > maxDispFactor * h) {
417 // std::cout << "Fallback to full fine tree for edge " << edge << std::endl;
418 EntityHandle feature_tri = 0;
419 CHKERR fineTree.closest_triangle(fineRoot, query, closest_point,
420 feature_tri);
421 if (feature_tri) {
422 disp = get_disp(closest_point);
423 }
424 }
425
426 if (disp > maxDispFactor * h) {
427 MOFEM_LOG("SELF", Sev::warning) <<
428 "The displacement for edge " << edge
429 << " is too large: " << disp << " (max allowed (-max_midnode_disp_factor * h) "
430 << maxDispFactor * h << "). Using original mid-edge position.";
431
432 } else {
433 midNodeCoord[0] = closest_point[0];
434 midNodeCoord[1] = closest_point[1];
435 midNodeCoord[2] = closest_point[2];
436 }
437
438 diffNodeCoord.resize(3);
439 auto t_diff_node_coord = getFTensor1FromPtr<3>(&diffNodeCoord(0));
440 t_diff_node_coord(i) = t_mid_node_coord(i) - t_ave_mid_coord(i);
441 dOf.resize(3);
442 auto t_dof = getFTensor1FromPtr<3>(&dOf(0));
443 t_dof(i) = t_diff_node_coord(i) / edge_shape_function_val;
444 if (entPtr->getNbOfCoeffs() > 3) {
445 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
446 "this method works only fields which are rank 3 or lower");
447 }
448 ent_data[0] = dOf[0];
449 ent_data[1] = dOf[1];
450 ent_data[2] = dOf[2];
451
453}
#define LOBATTO_PHI0(x)
Definitions taken from Hermes2d code.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ AINSWORTH_LOBATTO_BASE
Definition definitions.h:62
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
PetscErrorCode Legendre_polynomials(int p, double s, double *diff_s, double *L, double *diffL, int dim)
Calculate Legendre approximation basis.
#define MOFEM_LOG(channel, severity)
Log.
FTensor::Index< 'i', SPACE_DIM > i
double h
static constexpr double delta
virtual moab::Interface & get_moab()=0
boost::shared_ptr< FieldEntity > entPtr
Shared pointer to field entity data.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

◆ postProcess()

MoFEMErrorCode MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::postProcess ( )
virtual

Post-processing function executed at loop completion.

This virtual function is called once at the end of the loop. It is typically used for:

  • Assembling matrices and vectors to global system
  • Computing global variables (e.g., total internal energy)
  • Finalizing computation results
  • Cleanup operations

Example of iterating over DOFs:

for(_IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP_(this,"DISPLACEMENT",it)) {
// Process each DOF
}
Returns
Error code indicating success or failure

Reimplemented from MoFEM::BasicMethod.

Definition at line 455 of file Projection10NodeCoordsOnField.cpp.

◆ preProcess()

MoFEMErrorCode MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::preProcess ( )
virtual

Pre-processing function executed at loop initialization.

This virtual function is called once at the beginning of the loop. It is typically used for:

  • Zeroing matrices and vectors
  • Computing shape functions on reference elements
  • Preprocessing boundary conditions
  • Setting up initial computation state
Returns
Error code indicating success or failure

Reimplemented from MoFEM::BasicMethod.

Definition at line 209 of file Projection10NodeCoordsOnField.cpp.

209 {
211
212 // Build coarse skin triangles from volume mesh.
213 Range coarse_vols;
214 CHKERR mField.get_moab().get_entities_by_dimension(0, 3, coarse_vols);
215 if (coarse_vols.empty()) {
216 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
217 "no 3d entities found in coarse mesh");
218 }
219
220 Skinner coarse_skinner(&mField.get_moab());
221 coarseSkinTris.clear();
222 CHKERR coarse_skinner.find_skin(0, coarse_vols, false, coarseSkinTris);
223 if (coarseSkinTris.empty()) {
224 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
225 "no skin triangles found in coarse mesh");
226 }
227
228 coarseSkinEdges.clear();
229 CHKERR mField.get_moab().get_adjacencies(coarseSkinTris, 1, false,
231 moab::Interface::UNION);
232
233 if (fineTris.empty()) {
234 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
235 "no surface triangles found in fine mesh");
236 }
237
238 auto tri_normal = [&](const EntityHandle tri,
241 const EntityHandle *conn;
242 int num_nodes;
243 double coords[9];
244 CHKERR fineMoab.get_connectivity(tri, conn, num_nodes, true);
245 CHKERR fineMoab.get_coords(conn, num_nodes, coords);
246
247 CHKERR mField.getInterface<Tools>()->getTriNormal(coords, &normal(0));
248 normal.normalize();
250 };
251
252 const double fine_planar_cos =
253 std::cos(std::acos(planarAngleCos) * 1e-6);
254
255 Range fine_edges;
256 CHKERR fineMoab.get_adjacencies(fineTris, 1, true, fine_edges,
257 moab::Interface::UNION);
258
259 if (debugSurfaces) {
260 CHKERR CutMeshInterface::SaveData(fineMoab)("fine_edge_features.vtk",
261 fine_edges);
262 }
263 fineFeatureTris.clear();
264 for (auto edge : fine_edges) {
265 Range adj_tris;
266 CHKERR fineMoab.get_adjacencies(&edge, 1, 2, false, adj_tris,
267 moab::Interface::UNION);
268 adj_tris = intersect(adj_tris, fineTris);
269 if (adj_tris.size() != 2) {
270 continue;
271 }
272
275 CHKERR tri_normal(adj_tris[0], n0);
276 CHKERR tri_normal(adj_tris[1], n1);
277 const double dot = n0(i) * n1(i);
278 if (std::abs(dot) < fine_planar_cos) {
279 fineFeatureTris.insert(adj_tris[0]);
280 fineFeatureTris.insert(adj_tris[1]);
281 }
282 }
283
285 if (fineFeatureTris.empty()) {
287 }
288 if (debugSurfaces) {
289 CHKERR CutMeshInterface::SaveData(fineMoab)("fine_tri_features.vtk",
291 }
292
293 // KD tree on fine surface for closest-point queries.
294 fineRoot = 0;
295 CHKERR fineTree.build_tree(fineTris, &fineRoot);
296 fineFeatureRoot = 0;
298
300}
static MoFEMErrorCode filterOuterTrisByEdges(moab::Interface &moab, Range &tris, const int free_edge_threshold=1)

Member Data Documentation

◆ aveMidCoord

VectorDouble3 MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::aveMidCoord
protected

Definition at line 89 of file Projection10NodeCoordsOnField.hpp.

◆ coarseSkinEdges

Range MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::coarseSkinEdges
protected

Definition at line 86 of file Projection10NodeCoordsOnField.hpp.

◆ coarseSkinTris

Range MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::coarseSkinTris
protected

Definition at line 85 of file Projection10NodeCoordsOnField.hpp.

◆ coords

VectorDouble MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::coords
protected

Definition at line 88 of file Projection10NodeCoordsOnField.hpp.

◆ debugSurfaces

bool MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::debugSurfaces
protected

Definition at line 84 of file Projection10NodeCoordsOnField.hpp.

◆ diffNodeCoord

VectorDouble3 MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::diffNodeCoord
protected

Definition at line 91 of file Projection10NodeCoordsOnField.hpp.

◆ dOf

VectorDouble3 MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::dOf
protected

Definition at line 92 of file Projection10NodeCoordsOnField.hpp.

◆ fieldName

std::string MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::fieldName
protected

Definition at line 74 of file Projection10NodeCoordsOnField.hpp.

◆ fineFeatureRoot

EntityHandle MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::fineFeatureRoot
protected

Definition at line 80 of file Projection10NodeCoordsOnField.hpp.

◆ fineFeatureTree

AdaptiveKDTree MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::fineFeatureTree
protected

Definition at line 79 of file Projection10NodeCoordsOnField.hpp.

◆ fineFeatureTris

Range MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::fineFeatureTris
protected

Definition at line 81 of file Projection10NodeCoordsOnField.hpp.

◆ fineMoab

moab::Interface& MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::fineMoab
protected

Definition at line 72 of file Projection10NodeCoordsOnField.hpp.

◆ fineRoot

EntityHandle MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::fineRoot
protected

Definition at line 78 of file Projection10NodeCoordsOnField.hpp.

◆ fineTree

AdaptiveKDTree MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::fineTree
protected

Definition at line 77 of file Projection10NodeCoordsOnField.hpp.

◆ fineTris

Range MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::fineTris
protected

Definition at line 73 of file Projection10NodeCoordsOnField.hpp.

◆ maxDispFactor

double MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::maxDispFactor
protected

Definition at line 83 of file Projection10NodeCoordsOnField.hpp.

◆ mField

Interface& MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::mField
protected

Definition at line 71 of file Projection10NodeCoordsOnField.hpp.

◆ midNodeCoord

VectorDouble3 MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::midNodeCoord
protected

Definition at line 90 of file Projection10NodeCoordsOnField.hpp.

◆ planarAngleCos

double MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::planarAngleCos
protected

Definition at line 82 of file Projection10NodeCoordsOnField.hpp.

◆ vErbose

int MoFEM::Projection10NodeCoordsFromFineSurfaceOnField::vErbose
protected

Definition at line 75 of file Projection10NodeCoordsOnField.hpp.


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