OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
particletopologydescription.h
Go to the documentation of this file.
1 /*
2  *
3  * ##### ##### ###### ###### ### ###
4  * ## ## ## ## ## ## ## ### ##
5  * ## ## ## ## #### #### ## # ##
6  * ## ## ## ## ## ## ## ##
7  * ## ## ## ## ## ## ## ##
8  * ##### ##### ## ###### ## ##
9  *
10  *
11  * OOFEM : Object Oriented Finite Element Code
12  *
13  * Copyright (C) 1993 - 2013 Borek Patzak
14  *
15  *
16  *
17  * Czech Technical University, Faculty of Civil Engineering,
18  * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19  *
20  * This program is free software; you can redistribute it and/or modify
21  * it under the terms of the GNU General Public License as published by
22  * the Free Software Foundation; either version 2 of the License, or
23  * (at your option) any later version.
24  *
25  * This program is distributed in the hope that it will be useful,
26  * but WITHOUT ANY WARRANTY; without even the implied warranty of
27  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
28  * GNU General Public License for more details.
29  *
30  * You should have received a copy of the GNU General Public License
31  * along with this program; if not, write to the Free Software
32  * Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
33  */
34 
35 #ifndef particletopologydescription_h
36 #define particletopologydescription_h
37 
38 #include "topologydescription.h"
39 #include "intarray.h"
40 #include "floatarray.h"
41 #include "floatmatrix.h"
42 
43 #include <string>
44 #include <list>
45 #include <vector>
46 #include <memory>
47 
50 
51 #define _IFT_Circle_center "center"
52 #define _IFT_Circle_radius "radius"
53 #define _IFT_Circle_start "start"
54 #define _IFT_Circle_end "end"
55 
56 
58 
59 #define _IFT_Line_start "start"
60 #define _IFT_Line_end "end"
61 
62 
64 
65 #define _IFT_Point_coords "coords"
66 
67 
69 
70 #define _IFT_Meshing_elementType "elementtype"
71 #define _IFT_Meshing_set "set"
72 
73 
75 
76 #define _IFT_ParticleTopologyDescription_Name "particletopology"
77 #define _IFT_ParticleTopologyDescription_nsd "nsd"
78 #define _IFT_ParticleTopologyDescription_baseResolution "baseresolution"
79 #define _IFT_ParticleTopologyDescription_tubeWidth "tubewidth"
80 #define _IFT_ParticleTopologyDescription_neighbors "neighbors"
81 #define _IFT_ParticleTopologyDescription_boundingBoxA "bboxa"
82 #define _IFT_ParticleTopologyDescription_boundingBoxB "bboxb"
83 #define _IFT_ParticleTopologyDescription_numberOfSegments "nsegments"
84 #define _IFT_ParticleTopologyDescription_regionOutside "regionoutside"
85 #define _IFT_ParticleTopologyDescription_regionInside "regioninside"
86 #define _IFT_ParticleTopologyDescription_identification "id"
87 
88 
89 
90 namespace oofem
91 {
92 template< class Point >class ParticleGrid;
93 class Domain;
94 struct Triangle_PSLG;
95 
97 struct ParticlePoint {
99  ParticlePoint(const FloatArray &foot, int id, const FloatArray &normal, double distance2) :
100  foot(foot), id(id), normal(normal), distance2(distance2), removal(false) {};
101  FloatArray foot;
102  // Auxiliary information
104  int id;
107  double distance2;
108  int node;
109  int c_node;
110  bool removal;
111 };
112 
130 {
131 protected:
135  bool writeVTK;
137  bool resampled;
139  double maxdisp2;
141  double tubeWidth;
143  int m;
145  std :: unique_ptr< ParticleGrid< ParticlePoint > > grid;
146 
148  std :: list< ParticlePoint >corners;
149 
156  IntArray regionOutside, regionInside;
157  FloatMatrix mergeID, controlID;
158 
162  std :: vector< std :: string >regionElementType;
164 
169  void resample();
170 
176  TopologyState checkOverlap();
177 
183  void removePoints(ParticleGrid< ParticlePoint > &g) const;
184 
192  bool findDisplacement(FloatArray &answer, int id, const FloatArray &footpoint, TimeStep *tStep) const;
193 
202  void collectNeighbors(std :: list< ParticlePoint * > &answer, const ParticlePoint *p, double dist = 0) const;
203 
207  static void getBoundingBox(FloatArray &x0, FloatArray &x1, const FloatArray &c, double width);
208 
212  void calculateShortestDistance(const ParticlePoint *p, std :: list< ParticlePoint * > &points, ParticleGrid< ParticlePoint > &grid) const;
213 
217  double shortestDistanceFromCurve(const FloatArray &a, double txi_min, double txi_max,
218  const FloatArray &n0, const FloatArray &y0, const FloatArray &p, FloatArray &foot, FloatArray &normal) const;
219 
228  void addLineSegment(int id, const FloatArray &p0, const FloatArray &p1, ParticleGrid< ParticlePoint > &grid) const;
239  void addCircleSegment(int id, const FloatArray &c, double r, double v0, double v1, ParticleGrid< ParticlePoint > &grid) const;
246  void addCorner(int id, const FloatArray &c, ParticleGrid< ParticlePoint > &grid);
247 
252  void generatePSLG(Triangle_PSLG &PSLG);
253 
254 public:
256  virtual ~ParticleTopologyDescription();
257  virtual bool instanciateYourself(DataReader &dr);
258 
259  virtual TopologyState updateYourself(TimeStep *tStep);
260 
272  virtual void generateMesh(std :: vector< FloatArray > &nodes, std :: vector< IntArray > &elements, std :: vector< IntArray > &segments,
273  std :: vector< IntArray > &n_markers, IntArray &e_markers, IntArray &s_markers, IntArray &e_egt, IntArray &s_egt);
274 
275  virtual void replaceFEMesh();
276 
277  virtual void doOutput(TimeStep *tStep);
278  virtual void writeDataToFile(const char *name) const;
279  virtual void writeVTKFile(const char *name) const;
280 
281  virtual const char *giveClassName() const { return "ParticleTopologyDescription"; }
282 };
283 } // end namespace oofem
284 #endif // particletopologydescription_h
Class and object Domain.
Definition: domain.h:115
bool writeVTK
Conditional for printing VTK output.
IntArray regionOutside
Mapping of regions from delimited by each id.
FloatArray corner
Corner particle (for open surfaces).
Class representing the abstraction for input data source.
Definition: datareader.h:50
Particle grid data structure for n-D grids.
Definition: particlegrid.h:57
Class implementing an array of integers.
Definition: intarray.h:61
FloatArray total_displacement
Total displacement since last resampling.
FloatArray foot
Closest coordinate on surface.
bool resampled
Denotes if the active grid is newly resampled.
std::unique_ptr< ParticleGrid< ParticlePoint > > grid
The grid of points, the actual topological information.
ParticlePoint(const FloatArray &foot, int id, const FloatArray &normal, double distance2)
Default point type for describing topology.
int c_node
Corner node number (for open surfaces).
Plane straight line graph used as input for meshing with triangle.
Class representing vector of real numbers.
Definition: floatarray.h:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
int m
Number of points to use for resampling.
int node
Node number (for meshing)
Abstract class for topology description.
double tubeWidth
Width of the tube around the interfaces.
virtual const char * giveClassName() const
Gives the name of the class.
std::list< ParticlePoint > corners
Corner nodes.
double maxdisp2
Maximum squared displacement of any particle.
FloatArray normal
Surface normal at foot point.
the oofem namespace is to define a context or scope in which all oofem names are defined.
bool useDisplacements
Determines if velocity or displacements dofs should be used to update geometry.
std::vector< std::string > regionElementType
Mapping from region to FE components.
TopologyState
Determines the state of the evolving topology.
Class representing solution step.
Definition: timestep.h:80
double distance2
Squared distance stored for efficiency.
A grid based particle method for describing topology.

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:30 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011