File: frameSolver.cpp

package info (click to toggle)
gmsh 4.7.1%2Bds1-5
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 95,484 kB
  • sloc: cpp: 566,747; ansic: 150,384; yacc: 7,198; python: 6,130; java: 3,486; lisp: 622; lex: 621; makefile: 613; perl: 571; sh: 439; xml: 415; javascript: 113; pascal: 35; modula3: 32
file content (301 lines) | stat: -rw-r--r-- 9,465 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
// Gmsh - Copyright (C) 1997-2020 C. Geuzaine, J.-F. Remacle
//
// See the LICENSE.txt file for license information. Please report all
// issues on https://gitlab.onelab.info/gmsh/gmsh/issues.

#include "GmshConfig.h"
#include "GModel.h"
#include "GVertex.h"
#include "GEdge.h"
#include "frameSolver.h"
#include "linearSystemCSR.h"
#include "linearSystemPETSc.h"
#include "linearSystemFull.h"

#if defined(HAVE_POST)
#include "PView.h"
#include "PViewData.h"
#endif

frameSolver2d::frameSolver2d(GModel *gm) : pAssembler(0), _myModel(gm) {}

void frameSolver2d::addFixations(const std::vector<int> &dirs,
                                 const std::vector<int> &modelVertices,
                                 double value)
{
  for(std::size_t j = 0; j < modelVertices.size(); j++) {
    GVertex *gv = _myModel->getVertexByTag(modelVertices[j]);
    if(gv) {
      for(std::size_t i = 0; i < dirs.size(); i++) {
        _fixations.push_back(gmshFixation(gv, dirs[i], value));
      }
    }
  }
}

void frameSolver2d::addNodalForces(const std::vector<int> &modelVertices,
                                   const std::vector<double> &force)
{
  for(std::size_t j = 0; j < modelVertices.size(); j++) {
    GVertex *gv = _myModel->getVertexByTag(modelVertices[j]);
    if(gv) {
      _nodalForces.push_back(std::make_pair(gv, force));
    }
  }
}

void frameSolver2d::addBeamsOrBars(const std::vector<int> &modelEdges, double E,
                                   double I, double A, int r[2])
{
  int r_middle[2] = {1, 1}, r_left[2] = {r[0], 1}, r_right[2] = {0, r[1]};
  //  printf("adding %d beams\n",modelEdges.size());
  for(std::size_t i = 0; i < modelEdges.size(); i++) {
    GEdge *ge = _myModel->getEdgeByTag(modelEdges[i]);
    if(ge) {
      //      printf("model edge %d found\n",ge->tag());
      for(std::size_t j = 0; j < ge->lines.size(); ++j) {
        MLine *l = ge->lines[j];
        if(j == 0 && j == ge->lines.size() - 1)
          _beams.push_back(gmshBeam2d(l, E, I, A, r));
        else if(j == 0)
          _beams.push_back(gmshBeam2d(l, E, I, A, r_left));
        else if(j == ge->lines.size() - 1)
          _beams.push_back(gmshBeam2d(l, E, I, A, r_right));
        else
          _beams.push_back(gmshBeam2d(l, E, I, A, r_middle));
      }
    }
  }
}

void frameSolver2d::addBeams(const std::vector<int> &modelEdges, double E,
                             double I, double A)
{
  int r[2] = {1, 1};
  addBeamsOrBars(modelEdges, E, I, A, r);
}

void frameSolver2d::addBars(const std::vector<int> &modelEdges, double E,
                            double I, double A)
{
  int r[2] = {0, 0};
  addBeamsOrBars(modelEdges, E, I, A, r);
}

// solve any time a parameter is modified
/*

  A vertex is connected to beams. The question
  is how many bars are rotuled

  We define 2 dofs per node

 */
void frameSolver2d::createDofs()
{
  //  printf("coucou %d fixations\n",_fixations.size());
  for(std::size_t i = 0; i < _fixations.size(); ++i) {
    const gmshFixation &f = _fixations[i];
    //    printf("f._vertex(%d) = %p %d
    //    %g\n",i,f._vertex,f._direction,f._value);
    MVertex *v = f._vertex->mesh_vertices[0];
    Dof DOF(v->getNum(), f._direction);
    pAssembler->fixDof(DOF, f._value);
  }

  //  printf("coucou2\n");
  computeRotationTags();
  //  printf("coucou3\n");
  for(std::size_t i = 0; i < _beams.size(); i++) {
    //    printf("beam[%d] rot %d
    //    %d\n",i,_beams[i]._rotationTags[0],_beams[i]._rotationTags[1]);
    for(std::size_t j = 0; j < 2; j++) {
      MVertex *v = _beams[i]._element->getVertex(j);
      Dof theta(v->getNum(),
                Dof::createTypeWithTwoInts(2, _beams[i]._rotationTags[j]));
      pAssembler->numberDof(theta);
      Dof U(v->getNum(), 0);
      pAssembler->numberDof(U);
      Dof V(v->getNum(), 1);
      pAssembler->numberDof(V);
    }
  }
  //  printf("%d dofs\n",pAssembler->sizeOfR());
}

void frameSolver2d::computeStiffnessMatrix(int iBeam, fullMatrix<double> &K)
{
  const gmshBeam2d &b = _beams[iBeam];
  const double BS = b._e * b._i / (b._l * b._l * b._l);
  const double TS = b._e * b._a / b._l;
  const MVertex *v1 = b._element->getVertex(0);
  const MVertex *v2 = b._element->getVertex(1);
  const double alpha = atan2(v2->y() - v1->y(), v2->x() - v1->x());
  const double C = cos(alpha);
  const double S = sin(alpha);

  printf("beam %d %g %g %g\n", iBeam, alpha, C, S);

  fullMatrix<double> R(6, 6);
  R(0, 0) = R(1, 1) = R(3, 3) = R(4, 4) = C;
  R(0, 1) = R(3, 4) = S;
  R(1, 0) = R(4, 3) = -S;
  R(2, 2) = R(5, 5) = 1.0;

  fullMatrix<double> k(6, 6);

  // tensile stiffness
  k(0, 0) = k(3, 3) = TS;
  k(0, 3) = k(3, 0) = -TS;
  // bending stiffness
  k(1, 1) = k(4, 4) = 12 * BS;
  k(2, 2) = k(5, 5) = 4. * BS * b._l * b._l;
  k(1, 2) = k(2, 1) = k(1, 5) = k(5, 1) = 6 * BS * b._l;
  k(4, 2) = k(2, 4) = k(4, 5) = k(5, 4) = -6 * BS * b._l;
  k(4, 1) = k(1, 4) = -12 * BS;
  k(5, 2) = k(2, 5) = -2 * BS * b._l * b._l;

  fullMatrix<double> Rt(R), temp(6, 6);
  Rt.transposeInPlace();
  Rt.mult(k, temp);
  temp.mult(R, K);
}

void frameSolver2d::solve()
{
#if defined(HAVE_PETSC)
  linearSystemPETSc<double> *lsys = new linearSystemPETSc<double>;
#elif defined(HAVE_GMM)
  linearSystemCSRGmm<double> *lsys = new linearSystemCSRGmm<double>;
  lsys->setGmres(1);
  lsys->setNoisy(1);
#else
  linearSystemFull<double> *lsys = new linearSystemFull<double>;
#endif

  if(pAssembler) delete pAssembler;
  pAssembler = new dofManager<double>(lsys);

  // fix dofs and create free ones
  createDofs();

  // force vector
  std::vector<std::pair<GVertex *, std::vector<double> > >::iterator it =
    _nodalForces.begin();
  for(; it != _nodalForces.end(); ++it) {
    MVertex *v = it->first->mesh_vertices[0];
    const std::vector<double> &F = it->second;
    Dof DOFX(v->getNum(), 0);
    Dof DOFY(v->getNum(), 1);
    pAssembler->assemble(DOFX, F[0]);
    pAssembler->assemble(DOFY, F[1]);
  }

  // stifness matrix
  for(std::size_t i = 0; i < _beams.size(); i++) {
    fullMatrix<double> K(6, 6);
    computeStiffnessMatrix(i, K);
    _beams[i]._stiffness = K;
    MVertex *v0 = _beams[i]._element->getVertex(0);
    MVertex *v1 = _beams[i]._element->getVertex(1);
    Dof theta0(v0->getNum(),
               Dof::createTypeWithTwoInts(2, _beams[i]._rotationTags[0]));
    Dof theta1(v1->getNum(),
               Dof::createTypeWithTwoInts(2, _beams[i]._rotationTags[1]));
    Dof U0(v0->getNum(), 0);
    Dof U1(v1->getNum(), 0);
    Dof V0(v0->getNum(), 1);
    Dof V1(v1->getNum(), 1);
    Dof DOFS[6] = {U0, V0, theta0, U1, V1, theta1};
    for(int j = 0; j < 6; j++) {
      for(int k = 0; k < 6; k++) {
        pAssembler->assemble(DOFS[j], DOFS[k], K(j, k));
      }
    }
  }
  lsys->systemSolve();

  // save the solution
  for(std::size_t i = 0; i < _beams.size(); i++) {
    MVertex *v0 = _beams[i]._element->getVertex(0);
    MVertex *v1 = _beams[i]._element->getVertex(1);
    Dof theta0(v0->getNum(),
               Dof::createTypeWithTwoInts(2, _beams[i]._rotationTags[0]));
    Dof theta1(v1->getNum(),
               Dof::createTypeWithTwoInts(2, _beams[i]._rotationTags[1]));
    Dof U0(v0->getNum(), 0);
    Dof U1(v1->getNum(), 0);
    Dof V0(v0->getNum(), 1);
    Dof V1(v1->getNum(), 1);
    Dof DOFS[6] = {U0, V0, theta0, U1, V1, theta1};
    for(int j = 0; j < 6; j++) {
      pAssembler->getDofValue(DOFS[j], _beams[i]._displacement[j]);
    }
  }
  delete lsys;
  delete pAssembler;
}

void frameSolver2d::exportFrameData(const char *DISPL, const char *M)
{
#if defined(HAVE_POST)
  {
    std::map<int, std::vector<double> > data;
    for(std::size_t i = 0; i < _beams.size(); i++) {
      std::vector<double> tmp;
      // tmp.push_back(_beams[i]._e);
      // tmp.push_back(_beams[i]._i);
      // tmp.push_back(_beams[i]._a);
      tmp.reserve(6);
      for(int j = 0; j < 6; j++) {
        tmp.push_back(_beams[i]._displacement[j]);
      }
      data[_beams[i]._element->getNum()] = tmp;
    }
    PView *pv = new PView("displacements", "Beam", _myModel, data, 0.0, 6);
    pv->getData()->writeMSH(DISPL);
    delete pv;
  }
  {
    std::map<int, std::vector<double> > data;
    for(std::size_t i = 0; i < _beams.size(); i++) {
      std::vector<double> tmp;
      fullVector<double> d(_beams[i]._displacement, 6), F(6);
      _beams[i]._stiffness.mult(d, F);
      tmp.push_back(-F(2));
      tmp.push_back(F(5));
      data[_beams[i]._element->getNum()] = tmp;
    }
    PView *pv =
      new PView("Momentum", "ElementNodeData", _myModel, data, 0.0, 1);
    pv->getData()->writeMSH(M);
    delete pv;
  }
#endif
}

void frameSolver2d::computeRotationTags()
{
  std::multimap<MVertex *, gmshBeam2d *> v2b;
  for(std::size_t i = 0; i < _beams.size(); i++) {
    v2b.insert(std::make_pair(_beams[i]._element->getVertex(0), &_beams[i]));
    v2b.insert(std::make_pair(_beams[i]._element->getVertex(1), &_beams[i]));
  }

  std::multimap<MVertex *, gmshBeam2d *>::iterator s_it;
  for(std::multimap<MVertex *, gmshBeam2d *>::iterator it = v2b.begin();
      it != v2b.end(); it = s_it) {
    MVertex *theKey = it->first;

    std::pair<std::multimap<MVertex *, gmshBeam2d *>::iterator,
              std::multimap<MVertex *, gmshBeam2d *>::iterator>
      keyRange = v2b.equal_range(theKey);
    int countRotules = 0;
    for(s_it = keyRange.first; s_it != keyRange.second; ++s_it) {
      gmshBeam2d *b = s_it->second;
      if(!b->isRigid(theKey)) {
        b->setRotationTag(theKey, ++countRotules);
      }
    }
  }
}