-
Notifications
You must be signed in to change notification settings - Fork 0
/
small_strain_mech.cpp
1474 lines (1218 loc) · 42.7 KB
/
small_strain_mech.cpp
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
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/*
* Copyright (c) 2012-2015: G-CSC, Goethe University Frankfurt
* Authors: Raphael Prohl, Andreas Vogel
*
* This file is part of UG4.
*
* UG4 is free software: you can redistribute it and/or modify it under the
* terms of the GNU Lesser General Public License version 3 (as published by the
* Free Software Foundation) with the following additional attribution
* requirements (according to LGPL/GPL v3 §7):
*
* (1) The following notice must be displayed in the Appropriate Legal Notices
* of covered and combined works: "Based on UG4 (www.ug4.org/license)".
*
* (2) The following notice must be displayed at a prominent place in the
* terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
*
* (3) The following bibliography is recommended for citation and must be
* preserved in all covered files:
* "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
* parallel geometric multigrid solver on hierarchically distributed grids.
* Computing and visualization in science 16, 4 (2013), 151-164"
* "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
* flexible software system for simulating pde based models on high performance
* computers. Computing and visualization in science 16, 4 (2013), 165-179"
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU Lesser General Public License for more details.
*/
#include "small_strain_mech.h"
#include "material_laws/scaled_hooke_law.h"
// for various user data
#include "bindings/lua/lua_user_data.h"
#include "lib_disc/spatial_disc/user_data/user_data.h"
#include "lib_disc/spatial_disc/user_data/const_user_data.h"
#include "lib_disc/spatial_disc/disc_util/geom_provider.h"
#include "lib_disc/spatial_disc/disc_util/fe_geom.h"
#define PROFILE_SMALL_STRAIN_MECH
#ifdef PROFILE_SMALL_STRAIN_MECH
#define SMALL_STRAIN_MECH_PROFILE_FUNC() PROFILE_FUNC_GROUP("Small Strain Mech")
#define SMALL_STRAIN_MECH_PROFILE_BEGIN(name) PROFILE_BEGIN_GROUP(name, "Small Strain Mech")
#define SMALL_STRAIN_MECH_PROFILE_END() PROFILE_END()
#else
#define SMALL_STRAIN_MECH_PROFILE_FUNC()
#define SMALL_STRAIN_MECH_PROFILE_BEGIN(name)
#define SMALL_STRAIN_MECH_PROFILE_END()
#endif
using namespace std;
namespace ug {
namespace SmallStrainMechanics{
//////// Volume forces (IMPORT)
template<typename TDomain>
void SmallStrainMechanicsElemDisc<TDomain>::
set_volume_forces(SmartPtr<CplUserData<MathVector<dim>, dim> > user)
{
m_imVolForce.set_data(user);
}
template<typename TDomain>
void SmallStrainMechanicsElemDisc<TDomain>::
set_volume_forces(number vel_x)
{
SmartPtr<ConstUserVector<dim> > force(new ConstUserVector<dim>());
force->set_all_entries(vel_x);
set_volume_forces(force);
}
template<typename TDomain>
void SmallStrainMechanicsElemDisc<TDomain>::
set_volume_forces(number vel_x, number vel_y)
{
UG_THROW("SmallStrainMechanicsElemDisc: Setting force vector of dimension 2"
" to a Discretization for world dim " << dim);
}
template<>
void SmallStrainMechanicsElemDisc<Domain2d>::
set_volume_forces(number vel_x, number vel_y)
{
SmartPtr<ConstUserVector<dim> > force(new ConstUserVector<dim>());
force->set_entry(0, vel_x);
force->set_entry(1, vel_y);
set_volume_forces(force);
}
template<typename TDomain>
void SmallStrainMechanicsElemDisc<TDomain>::
set_volume_forces(number vel_x, number vel_y, number vel_z)
{
UG_THROW("SmallStrainMechanicsElemDisc: Setting force vector of dimension 3"
" to a Discretization for world dim " << dim);
}
template<>
void SmallStrainMechanicsElemDisc<Domain3d>::
set_volume_forces(number vel_x, number vel_y, number vel_z)
{
SmartPtr<ConstUserVector<dim> > force(new ConstUserVector<dim>());
force->set_entry(0, vel_x);
force->set_entry(1, vel_y);
force->set_entry(2, vel_z);
set_volume_forces(force);
}
#ifdef UG_FOR_LUA
template<typename TDomain>
void SmallStrainMechanicsElemDisc<TDomain>::
set_volume_forces(const char* fctName)
{
set_volume_forces(LuaUserDataFactory<MathVector<dim>, dim>::create(fctName));
}
#endif
/////////////// PRESSURE (IMPORT)
template<typename TDomain>
void SmallStrainMechanicsElemDisc<TDomain>::
set_div_factor(SmartPtr<CplUserData<number, dim> > user)
{
m_imDivergence.set_data(user);
}
template<typename TDomain>
void SmallStrainMechanicsElemDisc<TDomain>::
set_div_factor(number val)
{
set_div_factor(make_sp(new ConstUserNumber<dim>(val)));
}
#ifdef UG_FOR_LUA
template<typename TDomain>
void SmallStrainMechanicsElemDisc<TDomain>::
set_div_factor(const char* fctName)
{
set_div_factor(LuaUserDataFactory<number,dim>::create(fctName));
}
#endif
//////// Volume forces (IMPORT)
template<typename TDomain>
void SmallStrainMechanicsElemDisc<TDomain>::
set_viscous_forces(SmartPtr<CplUserData<MathVector<dim>, dim> > user1, SmartPtr<CplUserData<MathVector<dim>, dim> > user2)
{
m_imViscousForces[0].set_data(user1);
m_imViscousForces[1].set_data(user2);
}
template<typename TDomain>
void SmallStrainMechanicsElemDisc<TDomain>::
set_compress_factor(number val)
{
m_imCompressIndex.set_data(make_sp(new ConstUserNumber<dim>(val)));
}
//////////////////////////////////
template<typename TDomain>
void SmallStrainMechanicsElemDisc<TDomain>::
update_geo_elem(TBaseElem* elem, DimFEGeometry<dim>& geo)
{
SmartPtr<TDomain> dom = this->domain();
typedef typename IElemDisc<TDomain>::domain_type::position_accessor_type
position_accessor_type;
const position_accessor_type& aaPos = dom->position_accessor();
// coord and vertex array
MathVector<dim> coCoord[domain_traits<dim>::MaxNumVerticesOfElem];
// get vertices and extract corner coordinates
const size_t numVertices = elem->num_vertices();
for (size_t i = 0; i < numVertices; ++i) {
coCoord[i] = aaPos[elem->vertex(i)];
};
// prepare geometry for type and order
try{
geo.update(elem, &(coCoord[0]), m_lfeID, m_quadOrder);
}
UG_CATCH_THROW("SmallStrainMechanics::update_geo_elem:"
" Cannot update Finite Element Geometry.");
}
template<typename TDomain>
void
SmallStrainMechanicsElemDisc<TDomain>::
init_state_variables(const size_t order)
{
SMALL_STRAIN_MECH_PROFILE_BEGIN(SmallStrainMechInit_state_variables);
m_order = order;
m_lfeID = LFEID(LFEID::LAGRANGE, dim, order);
UG_LOG("\n");
// set default quadrature order if not set by user
if (!m_bQuadOrderUserDef) {
m_quadOrder = 2 * order + 1;
UG_LOG("Default QuadOrder is set in 'init_state_variables';"
"QuadOrder:" << m_quadOrder << "\n");
}
UG_LOG("In init_state: Order: " << m_order << " QuadOrder: "
<< m_quadOrder << "\n");
// clears and then attaches the attachments for the internal variables
// to the grid
SmartPtr<TDomain> dom = this->domain();
if(dom.invalid())
UG_THROW("Domain not provided for SmallStrainMechanics::init_state_variables");
typename TDomain::grid_type& grid = *(dom->grid());
if(m_spMatLaw.invalid())
UG_THROW("Material Law not provided for SmallStrainMechanics::init_state_variables");
m_spMatLaw->clear_attachments(grid);
m_spMatLaw->attach_internal_vars(grid);
// here the attachment is attached to the whole grid/multigrid
// and it remains attached until the
// destructor of SmallStrainMechanicsElemDisc!
DimFEGeometry<dim> geo;
typedef typename TDomain::grid_type::template traits<TBaseElem>::iterator
ElemIter;
for (ElemIter iter = grid.template begin<TBaseElem> (); iter
!= grid.template end<TBaseElem> (); iter++)
{
TBaseElem* elem = *iter;
update_geo_elem(elem, geo);
m_spMatLaw->init_internal_vars(elem, geo.num_ip());
}//end (ElemIter)
SMALL_STRAIN_MECH_PROFILE_END();
}
template<typename TDomain>
template<typename TElem, typename TFEGeom>
void SmallStrainMechanicsElemDisc<TDomain>::
prep_timestep_elem(const number time, const LocalVector& u,
GridObject* elem, const MathVector<dim> vCornerCoords[])
{
if (m_bOutWriter)
m_spOutWriter->pre_timestep();
}
template<typename TDomain>
template<typename TElem, typename TFEGeom>
void SmallStrainMechanicsElemDisc<TDomain>::
prep_elem_loop(const ReferenceObjectID roid, const int si)
{
// all this will be performed outside of the loop over the elements.
// Therefore it is not time critical.
// check, if a material law is set
if (m_spMatLaw.invalid())
UG_THROW("No material law set in "
"SmallStrainMechanicsElemDisc::prep_elem_loop \n");
// TODO: this only needs to be checked once at the beginning!
// -> prepare_setting?!
if (!m_spMatLaw->is_initialized())
m_spMatLaw->init();
// pass the material law to the output writer
if (m_bOutWriter && (!m_bMatLawPassedToOutWriter)){
m_spOutWriter->material_law(m_spMatLaw);
m_spOutWriter->quad_order(m_quadOrder);
m_bMatLawPassedToOutWriter = true;
}
// request geometry
TFEGeom& geo = GeomProvider<TFEGeom>::get(m_lfeID, m_quadOrder);
// prepare geometry for type and order
try{
geo.update_local(roid, m_lfeID, m_quadOrder);
}UG_CATCH_THROW("SmallStrainMechanicsElemDisc::prep_elem_loop:"
" Cannot update Finite Element Geometry.");
// set local positions for rhs
static const int refDim = TElem::dim;
m_imVolForce.template set_local_ips<refDim>(geo.local_ips(), geo.num_ip(), true);
m_imDivergence.template set_local_ips<refDim>(geo.local_ips(), geo.num_ip(), false);
m_imCompressIndex.template set_local_ips<refDim>(geo.local_ips(), geo.num_ip(), false);
m_imViscousForces[0].template set_local_ips<refDim>(geo.local_ips(), geo.num_ip(), true);
m_imViscousForces[1].template set_local_ips<refDim>(geo.local_ips(), geo.num_ip(), true);
}
template<typename TDomain>
template<typename TElem, typename TFEGeom>
void SmallStrainMechanicsElemDisc<TDomain>::fsh_elem_loop()
{}
template<typename TDomain>
template<typename TElem, typename TFEGeom>
void SmallStrainMechanicsElemDisc<TDomain>::
prep_elem(const LocalVector& u, GridObject* elem, const ReferenceObjectID roid, const MathVector<dim> vCornerCoords[])
{
// request geometry
TFEGeom& geo = GeomProvider<TFEGeom>::get(m_lfeID, m_quadOrder);
try{
geo.update(elem, vCornerCoords, m_lfeID, m_quadOrder);
}
UG_CATCH_THROW("SmallStrainMechanics::prep_elem:"
" Cannot update Finite Element Geometry.");
// set global positions for rhs
m_imVolForce.set_global_ips(geo.global_ips(), geo.num_ip());
m_imViscousForces[0].set_global_ips(geo.global_ips(), geo.num_ip());
m_imViscousForces[1].set_global_ips(geo.global_ips(), geo.num_ip());
m_imDivergence.set_global_ips(geo.global_ips(), geo.num_ip());
m_imCompressIndex.set_global_ips(geo.global_ips(), geo.num_ip());
// set pointer to internal variables of elem
m_spMatLaw->internal_vars(static_cast<TElem*>(elem));
}
// assemble stiffness jacobian
template<typename TDomain>
template<typename TElem, typename TFEGeom>
void SmallStrainMechanicsElemDisc<TDomain>::
add_jac_A_elem(LocalMatrix& J, const LocalVector& u,
GridObject* elem, const MathVector<dim> vCornerCoords[])
{
SMALL_STRAIN_MECH_PROFILE_BEGIN(SmallStrainMechAddJacA);
// request geometry
const TFEGeom& geo = GeomProvider<TFEGeom>::get(m_lfeID, m_quadOrder);
if (m_spMatLaw->elastTensIsConstant())
{
// material law with constant elasticity tensors
m_spElastTensor = m_spMatLaw->elasticityTensor();
}
/////// for brittle ///////
number scale = 1.0;
IScaledHookeLaw<TDomain>* pScaledHooke = dynamic_cast<IScaledHookeLaw<TDomain>*>(m_spMatLaw.get());
if(pScaledHooke != NULL){
scale = pScaledHooke->scaling_on_curr_elem();
}
/////// for brittle (end) ///////
MathMatrix<dim, dim> GradU;
for (size_t ip = 0; ip < geo.num_ip(); ++ip)
{
if (!(m_spMatLaw->elastTensIsConstant()))
{
// material law with non constant elasticity tensors
m_spMatLaw->template DisplacementGradient<TFEGeom>(GradU, ip, geo, u);
m_spElastTensor = m_spMatLaw->elasticityTensor(ip, geo.global_ip(ip), GradU);
}
// A) Compute Du:C:Dv = Du:sigma = sigma:Dv
for (size_t a = 0; a < geo.num_sh(); ++a) { // loop shape functions
// UG_LOG("add_jac_A_elem: sh="<<a);
for (size_t i = 0; i < (size_t) TDomain::dim; ++i) // loop component
for (size_t b = 0; b < geo.num_sh(); ++b) // loop shape functions
for (size_t j = 0; j < (size_t) TDomain::dim; ++j) // loop component
{
number integrandC = 0.0;
// Du:C:Dv = Du:sigma = sigma:Dv
for (size_t K = 0; K < (size_t) dim; ++K)
for (size_t L = 0; L < (size_t) dim; ++L)
{
integrandC += scale
* geo.global_grad(ip, a)[K]
* (*m_spElastTensor)[i][K][j][L]
* geo.global_grad(ip, b)[L];
}
J(i, a, j, b) += integrandC * geo.weight(ip);
} //end (j)
}
} //end(ip)
if(m_imCompressIndex.data_given()) {
for(size_t ip = 0; ip < geo.num_ip(); ++ip)
{ // for all integration points
for(size_t i2 = 0; i2 < num_fct(); ++i2){
for(size_t sh2 = 0; sh2 < geo.num_sh(); ++sh2)
{
for(size_t sh = 0; sh < geo.num_sh(); ++sh)
{ // for all shape functions
for(size_t i = 0; i < num_fct(); ++i)
{ // for all components (i=1,2,3)
J(i,sh,i2,sh2) += geo.weight(ip) * m_imCompressIndex[ip] *
geo.global_grad(ip, sh2)[i2] * geo.global_grad(ip, sh)[i];
}
}
}}
}
}
}
// assemble mass jacobian
template<typename TDomain>
template<typename TElem, typename TFEGeom>
void SmallStrainMechanicsElemDisc<TDomain>::
add_jac_M_elem(LocalMatrix& J, const LocalVector& u,
GridObject* elem, const MathVector<dim> vCornerCoords[])
{
if(m_spMatLaw->needs_to_add_jac_m())
{
// request geometry
const TFEGeom& geo = GeomProvider<TFEGeom>::get(m_lfeID, m_quadOrder);
for(size_t ip = 0; ip < geo.num_ip(); ++ip)
for(size_t i = 0; i < geo.num_sh(); ++i)
for(size_t j = 0; j < geo.num_sh(); ++j)
{
// same value for all components
number value = geo.shape(ip, i) * geo.shape(ip, j) * geo.weight(ip);
for(size_t c = 0; c < (size_t)dim; ++c)
{
J(c, i, c, j) += value * m_massScale;
}
}
}
}
// assemble stiffness defect
template<typename TDomain>
template<typename TElem, typename TFEGeom>
void SmallStrainMechanicsElemDisc<TDomain>::
add_def_A_elem(LocalVector& d, const LocalVector& u,
GridObject* elem, const MathVector<dim> vCornerCoords[])
{
SMALL_STRAIN_MECH_PROFILE_BEGIN(SmallStrainMechAddDefA);
// request geometry
const TFEGeom& geo = GeomProvider<TFEGeom>::get(m_lfeID, m_quadOrder);
////////////// for brittle //////////////////
number scale = 1.0;
number* energy = NULL;
IScaledHookeLaw<TDomain>* pScaledHooke = dynamic_cast<IScaledHookeLaw<TDomain>*>(m_spMatLaw.get());
if(pScaledHooke != NULL)
{
scale = pScaledHooke->scaling_on_curr_elem();
energy = &pScaledHooke->energy_on_curr_elem();
(*energy) = 0.0;
}
MathMatrix<dim,dim> strainTens;
number volElem = 0.0;
////////////// for brittle //////////////////
// a) Cauchy tensor
MathMatrix<dim, dim> sigma, GradU;
for (size_t ip = 0; ip < geo.num_ip(); ++ip)
{
// compute cauchy-stress tensor sigma at a ip
m_spMatLaw->template DisplacementGradient<TFEGeom>(GradU, ip, geo, u);
m_spMatLaw->stressTensor(sigma, ip, geo.global_ip(ip), GradU);
////////////// for brittle //////////////////
m_spElastTensor = m_spMatLaw->elasticityTensor(ip, geo.global_ip(ip), GradU);
if(pScaledHooke)
pScaledHooke->strainTensor(strainTens, GradU);
////////////// for brittle //////////////////
for (size_t sh = 0; sh < geo.num_sh(); ++sh)
{// loop shape functions
//UG_LOG("add_def_A_elem: sh="<<sh);
for (size_t i = 0; i < num_fct(); ++i) // loop components
{
number innerForcesIP = 0.0;
// i-th comp. of INTERNAL FORCES at node a:
for (size_t j = 0; j < (size_t) dim; ++j)
innerForcesIP += sigma[i][j] * geo.global_grad(ip, sh)[j];
//if (sh>=dim) {UG_LOG("add_def_A_elem: innerForcesIP="<<innerForcesIP << std::endl);}
d(i, sh) += geo.weight(ip) * scale * innerForcesIP;
////////////// for brittle //////////////////
if(energy){
for (size_t j = 0; j < (size_t) dim; ++j)
for (size_t K = 0; K < (size_t) dim; ++K)
for (size_t L = 0; L < (size_t) dim; ++L)
(*energy) += 0.5 * strainTens(i, K) * scale * (*m_spElastTensor)[i][K][j][L] * strainTens(j, L) * geo.weight(ip);
volElem += geo.weight(ip);
}
////////////// for brittle //////////////////
} //end (i)
}
}//end (ip)
////////////// for brittle //////////////////
if(energy){
(*energy) = (*energy) / volElem;
pScaledHooke->post_process_energy_on_curr_elem();
}
////////////// for brittle //////////////////
if(m_imCompressIndex.data_given()) {
for(size_t ip = 0; ip < geo.num_ip(); ++ip)
{ // for all integration points
number divU = 0.0;
for(size_t i = 0; i < num_fct(); ++i)
for(size_t sh = 0; sh < geo.num_sh(); ++sh)
{
divU += u(i,sh) * geo.global_grad(ip, sh)[i];
}
for(size_t sh = 0; sh < geo.num_sh(); ++sh)
{ // for all shape functions
for(size_t i = 0; i < num_fct(); ++i)
{ // for all components (i=1,2,3)
d(i,sh) += geo.weight(ip) * m_imCompressIndex[ip] *
divU * geo.global_grad(ip, sh)[i];
}
}
}
}
// scalar contribution, e.g, p * sum dx_i Phi_sh,i
if(m_imDivergence.data_given()) {
for(size_t sh = 0; sh < geo.num_sh(); ++sh)
{ // for all shape functions
for(size_t i = 0; i < num_fct(); ++i)
{ // for all components (i=1,2,3)
for(size_t ip = 0; ip < geo.num_ip(); ++ip)
{ // for all integration points
number ipGradI = geo.global_grad(ip, sh)[i];
d(i,sh) -= geo.weight(ip)*ipGradI*m_imDivergence[ip];
}
}
}
} // pressure data
}
// assemble mass-defect
template<typename TDomain>
template<typename TElem, typename TFEGeom>
void SmallStrainMechanicsElemDisc<TDomain>::
add_def_M_elem(LocalVector& d, const LocalVector& u,
GridObject* elem, const MathVector<dim> vCornerCoords[])
{}
// assemble right-hand-side d(i,sh)
template<typename TDomain>
template<typename TElem, typename TFEGeom>
void SmallStrainMechanicsElemDisc<TDomain>::
add_rhs_elem(LocalVector& d, GridObject* elem, const MathVector<dim> vCornerCoords[])
{
// request geometry
const TFEGeom& geo = GeomProvider<TFEGeom>::get(m_lfeID, m_quadOrder);
// a) volume forces, e.g. $$ \vec F*Phi =(grad p)*Phi$$
if(m_imVolForce.data_given()) {
// loop ip
for(size_t ip = 0; ip < geo.num_ip(); ++ip)
{
// loop shape
for(size_t sh = 0; sh < geo.num_sh(); ++sh)
{
// loop component
for(size_t i = 0; i < num_fct(); ++i)
{
d(i,sh) += geo.weight(ip) * geo.shape(ip, sh) * m_imVolForce[ip][i];
}
}
}
}
/*
// b) pressure part: p(ip) div (V_sh,ip)
// if (m_imDivergence.data_given()) {
for (size_t ip = 0; ip < geo.num_ip(); ++ip)
{ // loop IPs
for (size_t sh = 0; sh < geo.num_sh(); ++sh)
{ // loop shape functions
for (size_t i = 0; i < num_fct(); ++i)
{ // loop components
number divIP = 0.0;
for (size_t j = 0; j < (size_t) dim; ++j)
{
divIP += geo.global_grad(ip, sh)[j];
}
d(i, sh) += geo.weight(ip)*divIP; //import(ip)
} // end(i)
}// end (sh)
} // end (ip)
*/
// c) viscous stresses:
// $$ v0^T (grad Phi + grad Phi^T) v1 = v0^T (grad Phi) v1 + ((grad Phi) v0)^T v1 $$
if (m_imViscousForces[0].data_given() &&
m_imViscousForces[1].data_given())
{
for (size_t ip = 0; ip < geo.num_ip(); ++ip)
{ // loop IPs
for (size_t sh = 0; sh < geo.num_sh(); ++sh)
{ // loop shape functions
number gradPhi_v0 = VecDot(geo.global_grad(ip, sh), m_imViscousForces[0][ip]);
number gradPhi_v1 = VecDot(geo.global_grad(ip, sh), m_imViscousForces[1][ip]);
/*for (size_t i = 0; i < num_fct(); ++i)
{
// (grad Phi) v0
gradPhi_v0 += geo.global_grad(ip, sh)[i]*m_imViscousForces[0][ip][i];
// (grad Phi) v1
gradPhi_v1 += geo.global_grad(ip, sh)[i]*m_imViscousForces[1][ip][i];
}
*/
for (size_t i = 0; i < num_fct(); ++i)
{
d(i, sh) += geo.weight(ip)*m_imViscousForces[0][ip][i]*gradPhi_v1;
d(i, sh) += geo.weight(ip)*m_imViscousForces[1][ip][i]*gradPhi_v0;
} // end(i)
}// end (sh)
} // end (ip)
}
}
/** computes the linearized defect w.r.t to the pressure:
$$ \frac{\partial d(u,I_1), ... d(u,I_n)}{\partial I_i} $$
(since d is a vector, we obtain a bunch of vectors)
*/
template<typename TDomain>
template <typename TElem, typename TFEGeom>
void SmallStrainMechanicsElemDisc<TDomain>::
lin_def_pressure(const LocalVector& u,
std::vector<std::vector<number> > vvvLinDef[],
const size_t nip)
{
// request geometry
//static const typename TGeomProvider::Type& geo = TGeomProvider::get();
static const TFEGeom& geo = GeomProvider<TFEGeom>::get(m_lfeID, m_quadOrder);
UG_ASSERT(nip == geo.num_ip(), "Number of integration points does not match!");
for(size_t sh = 0; sh < geo.num_sh(); ++sh)
{ // loop test spaces
for (size_t i = 0; i < num_fct(); ++i)
{ // loop component
// loop integration points
for(size_t ip = 0; ip < geo.num_ip(); ++ip)
{
vvvLinDef[ip][i][sh] = -geo.weight(ip)*geo.global_grad(ip, sh)[i];
}
}
}
}
/** computes the linearized defect w.r.t to the volume forces
$$ \frac{\partial d(u,I_1), ... d(u,I_n)}{\partial I_i} $$
(since d is a vector, we obtain a bunch of diagonal matrices)
*/
template<typename TDomain>
template <typename TElem, typename TFEGeom>
void SmallStrainMechanicsElemDisc<TDomain>::
lin_def_volume_forces(const LocalVector& u,
std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
const size_t nip)
{
// request geometry
static const TFEGeom& geo = GeomProvider<TFEGeom>::get(m_lfeID, m_quadOrder);
UG_ASSERT(nip == geo.num_ip(), "Number of integration points does not match!");
for(size_t ip = 0; ip < geo.num_ip(); ++ip)
{ // loop integration points
for (size_t i = 0; i < num_fct(); ++i)
{ // loop component
for(size_t sh = 0; sh < geo.num_sh(); ++sh)
{ // loop test spaces
// add to local defect
vvvLinDef[ip][i][sh] = 0.0;
(vvvLinDef[ip][i][sh])[i] = geo.weight(ip)*geo.shape(ip, sh);
/*VecScale(vvvLinDef[ip][i][sh], geo.global_grad(ip, sh),
geo.weight(ip) * shape_ui);*/
}
}
}
}
/** computes the linearized defect w.r.t to the volume forces
$$ \frac{\partial d(u,I_1), ... d(u,I_n)}{\partial I_i} $$
(since d is a vector, we obtain a bunch of diagonal matrices)
*/
template<typename TDomain>
template <typename TElem, typename TFEGeom>
void SmallStrainMechanicsElemDisc<TDomain>::
lin_def_viscous_forces0(const LocalVector& u,
std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
const size_t nip)
{
// request geometry
static const TFEGeom& geo = GeomProvider<TFEGeom>::get(m_lfeID, m_quadOrder);
UG_ASSERT(nip == geo.num_ip(), "Number of integration points does not match!");
// loop IPs
for (size_t ip = 0; ip < geo.num_ip(); ++ip)
{
// loop shape functions
for (size_t sh = 0; sh < geo.num_sh(); ++sh)
{
//number gradPhi_v0 = VecDot(geo.global_grad(ip, sh), m_imViscousForces[0][ip]);
number gradPhi_v1 = VecDot(geo.global_grad(ip, sh), m_imViscousForces[1][ip]);
// loop components
for (size_t i = 0; i < num_fct(); ++i)
{
VecScale(vvvLinDef[ip][i][sh], m_imViscousForces[1][ip], geo.weight(ip)*geo.global_grad(ip, sh)[i]);
(vvvLinDef[ip][i][sh])[i] += geo.weight(ip)*gradPhi_v1;
} // end(i)
}// end (sh)
} // end (ip)end (ip)
}
/** computes the linearized defect w.r.t to the volume forces
$$ \frac{\partial d(u,I_1), ... d(u,I_n)}{\partial I_i} $$
(since d is a vector, we obtain a bunch of diagonal matrices)
*/
template<typename TDomain>
template <typename TElem, typename TFEGeom>
void SmallStrainMechanicsElemDisc<TDomain>::
lin_def_viscous_forces1(const LocalVector& u,
std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
const size_t nip)
{
// request geometry
static const TFEGeom& geo = GeomProvider<TFEGeom>::get(m_lfeID, m_quadOrder);
UG_ASSERT(nip == geo.num_ip(), "Number of integration points does not match!");
// loop IPs
for (size_t ip = 0; ip < geo.num_ip(); ++ip)
{
// loop shape functions
for (size_t sh = 0; sh < geo.num_sh(); ++sh)
{
number gradPhi_v0 = VecDot(geo.global_grad(ip, sh), m_imViscousForces[0][ip]);
//number gradPhi_v1 = VecDot(geo.global_grad(ip, sh), m_imViscousForces[1][ip]);
// loop components
for (size_t i = 0; i < num_fct(); ++i)
{
VecScale(vvvLinDef[ip][i][sh], m_imViscousForces[0][ip], geo.global_grad(ip, sh)[i]);
(vvvLinDef[ip][i][sh])[i] += gradPhi_v0;
} // end(i)
}// end (sh)
} // end (ip)end (ip)
}
template<typename TDomain>
template<typename TElem, typename TFEGeom>
void SmallStrainMechanicsElemDisc<TDomain>::
fsh_timestep_elem(const number time, const LocalVector& u,
GridObject* elem, const MathVector<dim> vCornerCoords[])
{
SMALL_STRAIN_MECH_PROFILE_BEGIN(SmallStrainMechFshTimeStepElem);
// request geometry
TFEGeom& geo = GeomProvider<TFEGeom>::get(m_lfeID, m_quadOrder);
try{
geo.update(elem, vCornerCoords, m_lfeID, m_quadOrder);
}
UG_CATCH_THROW("SmallStrainMechanics::fsh_timestep_elem:"
" Cannot update Finite Element Geometry.");
// pointer to internal variable of current elem
m_spMatLaw->internal_vars(static_cast<TElem*>(elem));
// call OutputWriter
if (m_bOutWriter){
SmartPtr<TDomain> dom = this->domain();
m_spOutWriter->post_timestep_elem(time, dom, geo, static_cast<TElem*>(elem), u);
}
// UPDATE PLASTIC STRAINS 'eps_p' AND HARDENING VARIABLE 'alpha'
MathMatrix<dim, dim> GradU;
for (size_t ip = 0; ip < geo.num_ip(); ++ip)
{
// update the internal (plastic) variables strain_p_new and alpha
m_spMatLaw->template DisplacementGradient<TFEGeom>(GradU, ip, geo, u);
m_spMatLaw->update_internal_vars(ip, GradU);
}
}
template <typename TDomain>
void
SmallStrainMechanicsElemDisc<TDomain>::
set_assemble_funcs()
{
// set default quadrature order if not set by user
if (!m_bQuadOrderUserDef) {
m_quadOrder = 2 * m_order + 1;
}
// set all non-set orders
else
{
if (m_quadOrder < 0){
m_quadOrder = 2 * m_order + 1;
}
}
register_all_fe_funcs(m_order, m_quadOrder);
}
////////////////////////////////////////////////////////////////////////////////
// Constructor
////////////////////////////////////////////////////////////////////////////////
template <typename TDomain>
SmallStrainMechanicsElemDisc<TDomain>::
SmallStrainMechanicsElemDisc(const char* functions, const char* subsets) :
IElemDisc<TDomain> (functions, subsets),
m_spMatLaw(SPNULL), m_spElastTensor(SPNULL), m_spOutWriter(SPNULL),
m_bOutWriter(false), m_bMatLawPassedToOutWriter(false),
m_exDivergence(new DataExport<number, dim>(functions)),
m_exDisplacement(new DataExport<MathVector<dim>, dim>(functions)),
m_exStressTensor(new DataExport<MathMatrix<dim,dim>, dim>(functions))
{
// check number of functions
if (this->num_fct() != (size_t) dim)
UG_THROW("Wrong number of functions: The ElemDisc 'SmallStrainMechanics'"
" needs exactly "<<dim<<" symbolic function.");
// register imports
this->register_import(m_imDivergence);
this->register_import(m_imCompressIndex);
this->register_import(m_imVolForce);
this->register_import(m_imViscousForces[0]);
this->register_import(m_imViscousForces[1]);
m_imVolForce.set_rhs_part();
m_imViscousForces[0].set_rhs_part();
m_imViscousForces[1].set_rhs_part();
m_imVolForce.set_comp_lin_defect(false);
// set defaults
m_order = 1;
m_bQuadOrderUserDef = false;
m_quadOrder = -1;
m_massScale = 1.0;
// update assemble functions
set_assemble_funcs();
}
template <typename TDomain>
SmallStrainMechanicsElemDisc<TDomain>::
~SmallStrainMechanicsElemDisc()
{
SmartPtr<TDomain> dom = this->domain();
typename TDomain::grid_type& grid = *(dom->grid());
m_spMatLaw->clear_attachments(grid);
}
////////////////////////////////////////////////////////////////////////////////
// export functions
////////////////////////////////////////////////////////////////////////////////
/// computes the stresses
template <typename TDomain>
template <typename TElem, typename TFEGeom>
void SmallStrainMechanicsElemDisc<TDomain>::
ex_stress_fe(MathMatrix<dim,dim> vValue[],
const MathVector<dim> vGlobIP[],
number time, int si,
const LocalVector& u,
GridObject* elem,
const MathVector<dim> vCornerCoords[],
const MathVector<TFEGeom::dim> vLocIP[],
const size_t nip,
bool bDeriv,
std::vector<std::vector<MathMatrix<dim,dim> > > vvvDeriv[])
{
// request geometry
const TFEGeom& geo = GeomProvider<TFEGeom>::get(m_lfeID, m_quadOrder);
// reference element
typedef typename reference_element_traits<TElem>::reference_element_type
ref_elem_type;
// reference dimension
static const int refDim = reference_element_traits<TElem>::dim;
// reference object id
static const ReferenceObjectID roid = ref_elem_type::REFERENCE_OBJECT_ID;
UG_ASSERT(bDeriv==false ,"Huhh: Some one was lazy - please implement derivatives!");
// a) Cauchy tensor
MathMatrix<dim, dim> GradU;
// FE
if(vLocIP == geo.local_ips())
{
// Loop ip
for(size_t ip = 0; ip < geo.num_ip(); ++ip)
{
MathMatrix<dim, dim> &sigmaIP = vValue[ip];
// compute cauchy-stress tensor sigma at a ip
m_spMatLaw->template DisplacementGradient<TFEGeom>(GradU, ip, geo, u);
m_spMatLaw->stressTensor(sigmaIP, ip, geo.global_ip(ip), GradU);
}
} else {
// general case
try{
// request for trial space
const LocalShapeFunctionSet<refDim>& rTrialSpace = LocalFiniteElementProvider::get<refDim>(roid, m_lfeID);
// number of shape functions
const size_t numSH = rTrialSpace.num_sh();
// storage for gradient function at ip
std::vector<MathVector<refDim> > vLocGrad(numSH);
MathVector<refDim> locGradA, globGradA;
// Reference Mapping
MathMatrix<dim, refDim> JTInv;
ReferenceMapping<ref_elem_type, dim> mapping(vCornerCoords);
// loop requested ips
for(size_t ip = 0; ip < nip; ++ip)
{
// computing Cauchy stress sigma:
MathMatrix<dim, dim> &sigmaIP = vValue[ip];
sigmaIP = 0.0;
// evaluate shape gradients at ip
rTrialSpace.grads(vLocGrad, vLocIP[ip]);
mapping.jacobian_transposed_inverse(JTInv, vLocIP[ip]);
// loop components
for (size_t a = 0; a < (size_t) dim; ++a)
{