Newer
Older
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
#ifdef HAVE_CONFIG_H
#include "config.h"
#endif
#include <dune/fufem/geometry/convexpolyhedron.hh>
#include <dune/contact/projections/normalprojection.hh>
#include "diameter.hh"
#include "multi-body-problem-data/bc.hh"
#include "multi-body-problem-data/cuboidgeometry.hh"
#include "multi-body-problem-data/myglobalfrictiondata.hh"
#include "multi-body-problem-data/weakpatch.hh"
#include "stackedblocksfactory.hh"
template <class GridType, int dims>
void StackedBlocksFactory<GridType, dims>::setBodies() {
// set up cuboid geometries
double const length = 1.00;
double const width = 0.27;
double const weakLength = 0.20;
std::vector<LocalVector> origins(this->bodyCount_);
std::vector<LocalVector> lowerWeakOrigins(this->bodyCount_);
std::vector<LocalVector> upperWeakOrigins(this->bodyCount_);
#if MY_DIM == 3
double const depth = 0.60;
origins[0] = {0,0,0};
lowerWeakOrigins[0] = {0.2,0,0};
upperWeakOrigins[0] = {0.2,width,0};
cuboidGeometries_[0] = std::make_shared<CuboidGeometry>(origins[0], lowerWeakOrigins[0], upperWeakOrigins[0], 0.0, weakLength, length, width, depth);
for (size_t i=1; i<this->bodyCount_-1; i++) {
origins[i] = cuboidGeometries_[i-1]->D;
lowerWeakOrigins[i] = {0.6,i*width,0};
upperWeakOrigins[i] = {0.6,(i+1)*width,0};
cuboidGeometries_[i] = std::make_shared<CuboidGeometry>(origins[i], lowerWeakOrigins[i], upperWeakOrigins[i], weakLength, weakLength, length, width, depth);
}
const size_t idx = this->bodyCount_-1;
origins[idx] = cuboidGeometries_[idx-1]->D;
lowerWeakOrigins[idx] = {0.6,idx*width,0};
upperWeakOrigins[idx] = {0.6,(idx+1)*width,0};
cuboidGeometries_[idx] = std::make_shared<CuboidGeometry>(origins[idx], lowerWeakOrigins[idx], upperWeakOrigins[idx], weakLength, 0.0, length, width, depth);
#elif MY_DIM == 2
origins[0] = {0,0};
lowerWeakOrigins[0] = {0.2,0};
upperWeakOrigins[0] = {0.2,width};
cuboidGeometries_[0] = std::make_shared<CuboidGeometry>(origins[0], lowerWeakOrigins[0], upperWeakOrigins[0], 0.0, weakLength, length, width);
for (size_t i=1; i<this->bodyCount_-1; i++) {
origins[i] = cuboidGeometries_[i-1]->D;
lowerWeakOrigins[i] = {0.6,i*width};
upperWeakOrigins[i] = {0.6,(i+1)*width};
cuboidGeometries_[i] = std::make_shared<CuboidGeometry>(origins[i], lowerWeakOrigins[i], upperWeakOrigins[i], weakLength, weakLength, length, width);
}
const size_t idx = this->bodyCount_-1;
origins[idx] = cuboidGeometries_[idx-1]->D;
lowerWeakOrigins[idx] = {0.6,idx*width};
upperWeakOrigins[idx] = {0.6,(idx+1)*width};
cuboidGeometries_[idx] = std::make_shared<CuboidGeometry>(origins[idx], lowerWeakOrigins[idx], upperWeakOrigins[idx], weakLength, 0.0, length, width);
#else
#error CuboidGeometry only supports 2D and 3D!"
#endif
// set up reference grids
gridConstructor_ = new GridsConstructor<GridType>(cuboidGeometries_);
auto& grids = gridConstructor_->getGrids();
for (size_t i=0; i<this->bodyCount_; i++) {
const auto& cuboidGeometry = *cuboidGeometries_[i];
// define weak patch and refine grid
auto lowerWeakPatch = lowerWeakPatches_[i] = std::make_shared<ConvexPolyhedron<LocalVector>>();
auto upperWeakPatch = upperWeakPatches_[i] = std::make_shared<ConvexPolyhedron<LocalVector>>();
getWeakPatch<LocalVector>(this->parset_.sub("boundary.friction.weakening"), cuboidGeometry, *lowerWeakPatch, *upperWeakPatch);
refine(*grids[i], *lowerWeakPatch, this->parset_.template get<double>("boundary.friction.smallestDiameter"), cuboidGeometry.lengthScale);
refine(*grids[i], *upperWeakPatch, this->parset_.template get<double>("boundary.friction.smallestDiameter"), cuboidGeometry.lengthScale);
// determine minDiameter and maxDiameter
double minDiameter = std::numeric_limits<double>::infinity();
double maxDiameter = 0.0;
for (auto &&e : elements(grids[i]->leafGridView())) {
auto const geometry = e.geometry();
auto const diam = diameter(geometry);
minDiameter = std::min(minDiameter, diam);
maxDiameter = std::max(maxDiameter, diam);
}
std::cout << "Grid" << i << " min diameter: " << minDiameter << std::endl;
std::cout << "Grid" << i << " max diameter: " << maxDiameter << std::endl;
}
for (size_t i=0; i<this->bodyCount_; i++) {
this->bodies_[i] = std::make_shared<typename Base::Body>(bodyData_, grids[i]);
}
}
template <class GridType, int dims>
void StackedBlocksFactory<GridType, dims>::setCouplings() {
const auto deformedGrids = this->levelContactNetwork_.deformedGrids();
for (size_t i=0; i<this->bodyCount_; i++) {
const auto& cuboidGeometry = *cuboidGeometries_[i];
leafFaces_[i] = std::make_shared<LeafFaces>(this->levelContactNetwork_.leafView(i), cuboidGeometry);
levelFaces_[i] = std::make_shared<LevelFaces>(this->levelContactNetwork_.levelView(i), cuboidGeometry);
}
auto contactProjection = std::make_shared<Dune::Contact::NormalProjection<LeafBoundaryPatch>>();
std::shared_ptr<typename Base::CouplingPair::BackEndType> backend = nullptr;
for (size_t i=0; i<this->couplings_.size(); i++) {
auto& coupling = this->couplings_[i];
coupling = std::make_shared<typename Base::CouplingPair>();
nonmortarPatch_[i] = std::make_shared<LevelBoundaryPatch>(levelFaces_[i]->upper);
mortarPatch_[i] = std::make_shared<LevelBoundaryPatch>(levelFaces_[i+1]->lower);
coupling->set(i, i+1, nonmortarPatch_[i], mortarPatch_[i], 0.1, Base::CouplingPair::CouplingType::STICK_SLIP, contactProjection, backend);
}
}
template <class GridType, int dims>
void StackedBlocksFactory<GridType, dims>::setBoundaryConditions() {
using LeafBoundaryCondition = BoundaryCondition<DeformedLeafGridView, dims>;
using FrictionBoundaryCondition = FrictionBoundaryCondition<DeformedLeafGridView, LocalVector, dims>;
using Function = Dune::VirtualFunction<double, double>;
std::shared_ptr<Function> neumannFunction = std::make_shared<NeumannCondition>();
std::shared_ptr<Function> velocityDirichletFunction = std::make_shared<VelocityDirichletCondition>();
const double lengthScale = CuboidGeometry::lengthScale;
for (size_t i=0; i<this->bodyCount_; i++) {
const auto& body = this->levelContactNetwork_.body(i);
const auto& leafVertexCount = body->leafVertexCount();
std::cout << "Grid" << i << " Number of DOFs: " << leafVertexCount << std::endl;
// friction boundary
if (i<this->bodyCount_-1 && upperWeakPatches_[i]->vertices.size()>0) {
std::shared_ptr<FrictionBoundaryCondition> frictionBoundary = std::make_shared<FrictionBoundaryCondition>();
frictionBoundary->setBoundaryPatch(leafFaces_[i]->upper);
frictionBoundary->setWeakeningPatch(upperWeakPatches_[i]);
frictionBoundary->setFrictionData(std::make_shared<MyGlobalFrictionData<LocalVector>>(this->parset_.sub("boundary.friction"), *upperWeakPatches_[i], lengthScale));
body->addBoundaryCondition(frictionBoundary);
}
if (i>0 && lowerWeakPatches_[i]->vertices.size()>0) {
std::shared_ptr<FrictionBoundaryCondition> frictionBoundary = std::make_shared<FrictionBoundaryCondition>();
frictionBoundary->setBoundaryPatch(leafFaces_[i]->lower);
frictionBoundary->setWeakeningPatch(lowerWeakPatches_[i]);
frictionBoundary->setFrictionData(std::make_shared<MyGlobalFrictionData<LocalVector>>(this->parset_.sub("boundary.friction"), *lowerWeakPatches_[i], lengthScale));
body->addBoundaryCondition(frictionBoundary);
}
// Neumann boundary
std::shared_ptr<LeafBoundaryCondition> neumannBoundary = std::make_shared<LeafBoundaryCondition>(std::make_shared<LeafBoundaryPatch>(body->leafView()), neumannFunction, "neumann");
body->addBoundaryCondition(neumannBoundary);
// Dirichlet Boundary
// identify Dirichlet nodes on leaf level
std::shared_ptr<Dune::BitSetVector<dims>> dirichletNodes = std::make_shared<Dune::BitSetVector<dims>>(leafVertexCount);
for (int j=0; j<leafVertexCount; j++) {
if (leafFaces_[i]->right.containsVertex(j))
(*dirichletNodes)[j][0] = true;
if (leafFaces_[i]->lower.containsVertex(j))
(*dirichletNodes)[j][0] = true;
#if MY_DIM == 3
if (leafFaces_[i]->front.containsVertex(j) || leafFaces_[i]->back.containsVertex(j))
dirichletNodes->at(j)[2] = true;
#endif
}
std::shared_ptr<LeafBoundaryCondition> dirichletBoundary = std::make_shared<LeafBoundaryCondition>("dirichlet");
dirichletBoundary->setBoundaryPatch(body->leafView(), dirichletNodes);
dirichletBoundary->setBoundaryFunction(velocityDirichletFunction);
body->addBoundaryCondition(dirichletBoundary);
}
}
#include "stackedblocksfactory_tmpl.cc"