gmsh-TingyuanDoc  0.1
An Open-Source Timing-driven Analytical Mixed-size FPGA Placer
frameSolver.cpp
Go to the documentation of this file.
1 // Gmsh - Copyright (C) 1997-2022 C. Geuzaine, J.-F. Remacle
2 //
3 // See the LICENSE.txt file in the Gmsh root directory for license information.
4 // Please report all issues on https://gitlab.onelab.info/gmsh/gmsh/issues.
5 
6 #include "GmshConfig.h"
7 #include "GModel.h"
8 #include "GVertex.h"
9 #include "GEdge.h"
10 #include "frameSolver.h"
11 #include "linearSystemCSR.h"
12 #include "linearSystemPETSc.h"
13 #include "linearSystemFull.h"
14 
15 #if defined(HAVE_POST)
16 #include "PView.h"
17 #include "PViewData.h"
18 #endif
19 
20 frameSolver2d::frameSolver2d(GModel *gm) : pAssembler(nullptr), _myModel(gm) {}
21 
22 void frameSolver2d::addFixations(const std::vector<int> &dirs,
23  const std::vector<int> &modelVertices,
24  double value)
25 {
26  for(std::size_t j = 0; j < modelVertices.size(); j++) {
27  GVertex *gv = _myModel->getVertexByTag(modelVertices[j]);
28  if(gv) {
29  for(std::size_t i = 0; i < dirs.size(); i++) {
30  _fixations.push_back(gmshFixation(gv, dirs[i], value));
31  }
32  }
33  }
34 }
35 
36 void frameSolver2d::addNodalForces(const std::vector<int> &modelVertices,
37  const std::vector<double> &force)
38 {
39  for(std::size_t j = 0; j < modelVertices.size(); j++) {
40  GVertex *gv = _myModel->getVertexByTag(modelVertices[j]);
41  if(gv) {
42  _nodalForces.push_back(std::make_pair(gv, force));
43  }
44  }
45 }
46 
47 void frameSolver2d::addBeamsOrBars(const std::vector<int> &modelEdges, double E,
48  double I, double A, int r[2])
49 {
50  int r_middle[2] = {1, 1}, r_left[2] = {r[0], 1}, r_right[2] = {0, r[1]};
51  // printf("adding %d beams\n",modelEdges.size());
52  for(std::size_t i = 0; i < modelEdges.size(); i++) {
53  GEdge *ge = _myModel->getEdgeByTag(modelEdges[i]);
54  if(ge) {
55  // printf("model edge %d found\n",ge->tag());
56  for(std::size_t j = 0; j < ge->lines.size(); ++j) {
57  MLine *l = ge->lines[j];
58  if(j == 0 && j == ge->lines.size() - 1)
59  _beams.push_back(gmshBeam2d(l, E, I, A, r));
60  else if(j == 0)
61  _beams.push_back(gmshBeam2d(l, E, I, A, r_left));
62  else if(j == ge->lines.size() - 1)
63  _beams.push_back(gmshBeam2d(l, E, I, A, r_right));
64  else
65  _beams.push_back(gmshBeam2d(l, E, I, A, r_middle));
66  }
67  }
68  }
69 }
70 
71 void frameSolver2d::addBeams(const std::vector<int> &modelEdges, double E,
72  double I, double A)
73 {
74  int r[2] = {1, 1};
75  addBeamsOrBars(modelEdges, E, I, A, r);
76 }
77 
78 void frameSolver2d::addBars(const std::vector<int> &modelEdges, double E,
79  double I, double A)
80 {
81  int r[2] = {0, 0};
82  addBeamsOrBars(modelEdges, E, I, A, r);
83 }
84 
85 // solve any time a parameter is modified
86 /*
87 
88  A vertex is connected to beams. The question
89  is how many bars are rotuled
90 
91  We define 2 dofs per node
92 
93  */
95 {
96  // printf("coucou %d fixations\n",_fixations.size());
97  for(std::size_t i = 0; i < _fixations.size(); ++i) {
98  const gmshFixation &f = _fixations[i];
99  // printf("f._vertex(%d) = %p %d
100  // %g\n",i,f._vertex,f._direction,f._value);
101  MVertex *v = f._vertex->mesh_vertices[0];
102  Dof DOF(v->getNum(), f._direction);
103  pAssembler->fixDof(DOF, f._value);
104  }
105 
106  // printf("coucou2\n");
108  // printf("coucou3\n");
109  for(std::size_t i = 0; i < _beams.size(); i++) {
110  // printf("beam[%d] rot %d
111  // %d\n",i,_beams[i]._rotationTags[0],_beams[i]._rotationTags[1]);
112  for(std::size_t j = 0; j < 2; j++) {
113  MVertex *v = _beams[i]._element->getVertex(j);
114  Dof theta(v->getNum(),
115  Dof::createTypeWithTwoInts(2, _beams[i]._rotationTags[j]));
116  pAssembler->numberDof(theta);
117  Dof U(v->getNum(), 0);
118  pAssembler->numberDof(U);
119  Dof V(v->getNum(), 1);
120  pAssembler->numberDof(V);
121  }
122  }
123  // printf("%d dofs\n",pAssembler->sizeOfR());
124 }
125 
127 {
128  const gmshBeam2d &b = _beams[iBeam];
129  const double BS = b._e * b._i / (b._l * b._l * b._l);
130  const double TS = b._e * b._a / b._l;
131  const MVertex *v1 = b._element->getVertex(0);
132  const MVertex *v2 = b._element->getVertex(1);
133  const double alpha = atan2(v2->y() - v1->y(), v2->x() - v1->x());
134  const double C = cos(alpha);
135  const double S = sin(alpha);
136 
137  printf("beam %d %g %g %g\n", iBeam, alpha, C, S);
138 
139  fullMatrix<double> R(6, 6);
140  R(0, 0) = R(1, 1) = R(3, 3) = R(4, 4) = C;
141  R(0, 1) = R(3, 4) = S;
142  R(1, 0) = R(4, 3) = -S;
143  R(2, 2) = R(5, 5) = 1.0;
144 
145  fullMatrix<double> k(6, 6);
146 
147  // tensile stiffness
148  k(0, 0) = k(3, 3) = TS;
149  k(0, 3) = k(3, 0) = -TS;
150  // bending stiffness
151  k(1, 1) = k(4, 4) = 12 * BS;
152  k(2, 2) = k(5, 5) = 4. * BS * b._l * b._l;
153  k(1, 2) = k(2, 1) = k(1, 5) = k(5, 1) = 6 * BS * b._l;
154  k(4, 2) = k(2, 4) = k(4, 5) = k(5, 4) = -6 * BS * b._l;
155  k(4, 1) = k(1, 4) = -12 * BS;
156  k(5, 2) = k(2, 5) = -2 * BS * b._l * b._l;
157 
158  fullMatrix<double> Rt(R), temp(6, 6);
159  Rt.transposeInPlace();
160  Rt.mult(k, temp);
161  temp.mult(R, K);
162 }
163 
165 {
166 #if defined(HAVE_PETSC)
168 #elif defined(HAVE_GMM)
170  lsys->setGmres(1);
171  lsys->setNoisy(1);
172 #else
174 #endif
175 
176  if(pAssembler) delete pAssembler;
177  pAssembler = new dofManager<double>(lsys);
178 
179  // fix dofs and create free ones
180  createDofs();
181 
182  // force vector
183  auto it =
184  _nodalForces.begin();
185  for(; it != _nodalForces.end(); ++it) {
186  MVertex *v = it->first->mesh_vertices[0];
187  const std::vector<double> &F = it->second;
188  Dof DOFX(v->getNum(), 0);
189  Dof DOFY(v->getNum(), 1);
190  pAssembler->assemble(DOFX, F[0]);
191  pAssembler->assemble(DOFY, F[1]);
192  }
193 
194  // stifness matrix
195  for(std::size_t i = 0; i < _beams.size(); i++) {
196  fullMatrix<double> K(6, 6);
198  _beams[i]._stiffness = K;
199  MVertex *v0 = _beams[i]._element->getVertex(0);
200  MVertex *v1 = _beams[i]._element->getVertex(1);
201  Dof theta0(v0->getNum(),
202  Dof::createTypeWithTwoInts(2, _beams[i]._rotationTags[0]));
203  Dof theta1(v1->getNum(),
204  Dof::createTypeWithTwoInts(2, _beams[i]._rotationTags[1]));
205  Dof U0(v0->getNum(), 0);
206  Dof U1(v1->getNum(), 0);
207  Dof V0(v0->getNum(), 1);
208  Dof V1(v1->getNum(), 1);
209  Dof DOFS[6] = {U0, V0, theta0, U1, V1, theta1};
210  for(int j = 0; j < 6; j++) {
211  for(int k = 0; k < 6; k++) {
212  pAssembler->assemble(DOFS[j], DOFS[k], K(j, k));
213  }
214  }
215  }
216  lsys->systemSolve();
217 
218  // save the solution
219  for(std::size_t i = 0; i < _beams.size(); i++) {
220  MVertex *v0 = _beams[i]._element->getVertex(0);
221  MVertex *v1 = _beams[i]._element->getVertex(1);
222  Dof theta0(v0->getNum(),
223  Dof::createTypeWithTwoInts(2, _beams[i]._rotationTags[0]));
224  Dof theta1(v1->getNum(),
225  Dof::createTypeWithTwoInts(2, _beams[i]._rotationTags[1]));
226  Dof U0(v0->getNum(), 0);
227  Dof U1(v1->getNum(), 0);
228  Dof V0(v0->getNum(), 1);
229  Dof V1(v1->getNum(), 1);
230  Dof DOFS[6] = {U0, V0, theta0, U1, V1, theta1};
231  for(int j = 0; j < 6; j++) {
232  pAssembler->getDofValue(DOFS[j], _beams[i]._displacement[j]);
233  }
234  }
235  delete lsys;
236  delete pAssembler;
237 }
238 
239 void frameSolver2d::exportFrameData(const char *DISPL, const char *M)
240 {
241 #if defined(HAVE_POST)
242  {
243  std::map<int, std::vector<double> > data;
244  for(std::size_t i = 0; i < _beams.size(); i++) {
245  std::vector<double> tmp;
246  // tmp.push_back(_beams[i]._e);
247  // tmp.push_back(_beams[i]._i);
248  // tmp.push_back(_beams[i]._a);
249  tmp.reserve(6);
250  for(int j = 0; j < 6; j++) {
251  tmp.push_back(_beams[i]._displacement[j]);
252  }
253  data[_beams[i]._element->getNum()] = tmp;
254  }
255  PView *pv = new PView("displacements", "Beam", _myModel, data, 0.0, 6);
256  pv->getData()->writeMSH(DISPL);
257  delete pv;
258  }
259  {
260  std::map<int, std::vector<double> > data;
261  for(std::size_t i = 0; i < _beams.size(); i++) {
262  std::vector<double> tmp;
263  fullVector<double> d(_beams[i]._displacement, 6), F(6);
264  _beams[i]._stiffness.mult(d, F);
265  tmp.push_back(-F(2));
266  tmp.push_back(F(5));
267  data[_beams[i]._element->getNum()] = tmp;
268  }
269  PView *pv =
270  new PView("Momentum", "ElementNodeData", _myModel, data, 0.0, 1);
271  pv->getData()->writeMSH(M);
272  delete pv;
273  }
274 #endif
275 }
276 
278 {
279  std::multimap<MVertex *, gmshBeam2d *> v2b;
280  for(std::size_t i = 0; i < _beams.size(); i++) {
281  v2b.insert(std::make_pair(_beams[i]._element->getVertex(0), &_beams[i]));
282  v2b.insert(std::make_pair(_beams[i]._element->getVertex(1), &_beams[i]));
283  }
284 
285  std::multimap<MVertex *, gmshBeam2d *>::iterator s_it;
286  for(auto it = v2b.begin(); it != v2b.end(); it = s_it) {
287  MVertex *theKey = it->first;
288 
289  std::pair<std::multimap<MVertex *, gmshBeam2d *>::iterator,
290  std::multimap<MVertex *, gmshBeam2d *>::iterator>
291  keyRange = v2b.equal_range(theKey);
292  int countRotules = 0;
293  for(s_it = keyRange.first; s_it != keyRange.second; ++s_it) {
294  gmshBeam2d *b = s_it->second;
295  if(!b->isRigid(theKey)) {
296  b->setRotationTag(theKey, ++countRotules);
297  }
298  }
299  }
300 }
linearSystemFull
Definition: linearSystemFull.h:17
PView
Definition: PView.h:27
frameSolver.h
fullVector< double >
linearSystemPETSc
Definition: linearSystemPETSc.h:150
MLine::getVertex
virtual MVertex * getVertex(int num)
Definition: MLine.h:45
GEdge::lines
std::vector< MLine * > lines
Definition: GEdge.h:42
F
#define F
Definition: DefaultOptions.h:23
MVertex
Definition: MVertex.h:24
linearSystemPETSc.h
LegendrePolynomials::f
void f(int n, double u, double *val)
Definition: orthogonalBasis.cpp:77
MVertex::getNum
std::size_t getNum() const
Definition: MVertex.h:86
frameSolver2d::pAssembler
dofManager< double > * pAssembler
Definition: frameSolver.h:58
frameSolver2d::_myModel
GModel * _myModel
Definition: frameSolver.h:62
PView.h
linearSystemCSR.h
gmshBeam2d::_i
double _i
Definition: frameSolver.h:20
gmshBeam2d
Definition: frameSolver.h:18
frameSolver2d::_fixations
std::vector< gmshFixation > _fixations
Definition: frameSolver.h:61
frameSolver2d::addBars
void addBars(const std::vector< int > &modelEdges, double E, double I, double A)
Definition: frameSolver.cpp:78
GModel::getEdgeByTag
GEdge * getEdgeByTag(int n) const
Definition: GModel.cpp:336
PViewData.h
MLine
Definition: MLine.h:21
linearSystemCSRGmm::setGmres
void setGmres(int n)
Definition: linearSystemCSR.h:189
frameSolver2d::_nodalForces
std::vector< std::pair< GVertex *, std::vector< double > > > _nodalForces
Definition: frameSolver.h:60
fullMatrix< double >
GEdge.h
gmshBeam2d::_element
MLine * _element
Definition: frameSolver.h:19
Dof
Definition: dofManager.h:19
GVertex
Definition: GVertex.h:23
frameSolver2d::frameSolver2d
frameSolver2d(GModel *myModel)
Definition: frameSolver.cpp:20
Dof::createTypeWithTwoInts
static int createTypeWithTwoInts(int i1, int i2)
Definition: dofManager.h:28
linearSystemFull::systemSolve
virtual int systemSolve()
Definition: linearSystemFull.h:77
GModel
Definition: GModel.h:44
PView::getData
PViewData * getData(bool useAdaptiveIfAvailable=false)
Definition: PView.cpp:233
frameSolver2d::_beams
std::vector< gmshBeam2d > _beams
Definition: frameSolver.h:59
gmshBeam2d::isRigid
bool isRigid(MVertex *v) const
Definition: frameSolver.h:35
PViewData::writeMSH
virtual bool writeMSH(const std::string &fileName, double version=2.2, bool binary=false, bool saveMesh=true, bool multipleView=false, int partitionNum=-1, bool saveInterpolationMatrices=true, bool forceNodeData=false, bool forceElementData=false)
Definition: PViewDataIO.cpp:202
dofManager::numberDof
virtual void numberDof(Dof key)
Definition: dofManager.h:213
dofManager< double >
frameSolver2d::computeRotationTags
void computeRotationTags()
Definition: frameSolver.cpp:277
linearSystemFull.h
dofManager::assemble
virtual void assemble(const Dof &R, const Dof &C, const dataMat &value)
Definition: dofManager.h:393
frameSolver2d::exportFrameData
void exportFrameData(const char *displ, const char *M)
Definition: frameSolver.cpp:239
frameSolver2d::addBeamsOrBars
void addBeamsOrBars(const std::vector< int > &modelEdges, double E, double I, double A, int r[2])
Definition: frameSolver.cpp:47
linearSystemCSRGmm
Definition: linearSystemCSR.h:176
S
#define S
Definition: DefaultOptions.h:21
dofManager::getDofValue
virtual void getDofValue(std::vector< Dof > &keys, std::vector< dataVec > &Vals)
Definition: dofManager.h:235
frameSolver2d::addBeams
void addBeams(const std::vector< int > &modelEdges, double E, double I, double A)
Definition: frameSolver.cpp:71
gmshFixation
Definition: frameSolver.h:46
GVertex.h
GEdge
Definition: GEdge.h:26
frameSolver2d::addFixations
void addFixations(const std::vector< int > &dirs, const std::vector< int > &modelVertices, double value)
Definition: frameSolver.cpp:22
frameSolver2d::solve
void solve()
Definition: frameSolver.cpp:164
GModel.h
fullMatrix::mult
void mult(const fullVector< scalar > &x, fullVector< scalar > &y) const
Definition: fullMatrix.h:487
MVertex::y
double y() const
Definition: MVertex.h:61
dofManager::fixDof
virtual void fixDof(Dof key, const dataVec &value)
Definition: dofManager.h:153
frameSolver2d::addNodalForces
void addNodalForces(const std::vector< int > &modelVertices, const std::vector< double > &force)
Definition: frameSolver.cpp:36
linearSystemCSRGmm::setNoisy
void setNoisy(int n)
Definition: linearSystemCSR.h:188
gmshBeam2d::_a
double _a
Definition: frameSolver.h:20
GModel::getVertexByTag
GVertex * getVertexByTag(int n) const
Definition: GModel.cpp:346
gmshBeam2d::setRotationTag
void setRotationTag(MVertex *v, int tag)
Definition: frameSolver.h:39
gmshBeam2d::_e
double _e
Definition: frameSolver.h:20
frameSolver2d::createDofs
void createDofs()
Definition: frameSolver.cpp:94
MVertex::x
double x() const
Definition: MVertex.h:60
frameSolver2d::computeStiffnessMatrix
void computeStiffnessMatrix(int iBeam, fullMatrix< double > &K)
Definition: frameSolver.cpp:126
gmshBeam2d::_l
double _l
Definition: frameSolver.h:20