Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
LinElas_def.hpp
1#ifndef LINELAS_def_hpp
2#define LINELAS_def_hpp
3#include "LinElas_decl.hpp"
4namespace FEDD {
5// Funktion fuer den homogenen Dirichletrand
6// void ZeroDirichlet(double* x, double* result, double t, double* parameters)
7// {
8// result[0] = 0.;
9// return;
10// }
11
12// TODO: Fuer FSI auf 0 und 0 abaendern!!!!!!!!!!
13// Funktion fuer die rechte Seite der DGL in 2D
14void rhsFunc2D(double* x, double* result, double* parameters)
15{
16 // Wir setzen die rechte Seite g_vec als g_vec = (0, g), mit g = -2.
17 result[0] = 0.;
18 result[1] = 0.;
19// if (parameters[0]<0.1) {
20// result[1] = -.2;
21// }
22
23 return;
24}
25
26// Funktion fuer die rechte Seite der DGL in 3D
27void rhsFunc3D(double* x, double* result, double* parameters)
28{
29 // Wir setzen die rechte Seite g_vec als g_vec = (0, g, 0), mit g = -2.
30 result[0] = 0.;
31 // result[1] = -2.;
32 result[1] = 0.;
33 result[2] = 0.;
34 return;
35}
36
37
38template<class SC,class LO,class GO,class NO>
39LinElas<SC,LO,GO,NO>::LinElas(const DomainConstPtr_Type &domain, std::string FEType, ParameterListPtr_Type parameterList):
40Problem_Type(parameterList, domain->getComm() )
41{
42 // Bem.: Hier benutzen wir auch direkt den Konstruktor der Klasse Problem, wo z.B. parameterList auf parameterList_ gesetzt wird
43
44 // d_s steht displacement (d) der Struktur (_s). Im Allgemeinen vektorwertig, daher domainDisplacement->GetDimension().
45 // Andernfalls benutzen wir stattdessen 1
46
47 // Siehe Problem_def.hpp fuer die Funktion AddVariable()
48 this->addVariable( domain , FEType , "d_s" ,domain->getDimension() );
49 this->dim_ = this->getDomain(0)->getDimension();
50}
51
52template<class SC,class LO,class GO,class NO>
53LinElas<SC,LO,GO,NO>::~LinElas()
54{
55
56}
57
58template<class SC,class LO,class GO,class NO>
59void LinElas<SC,LO,GO,NO>::info(){
60 this->infoProblem();
61}
62
63template<class SC,class LO,class GO,class NO>
64void LinElas<SC,LO,GO,NO>::assemble( std::string type ) const
65{
66 if(this->verbose_)
67 std::cout << "-- Assembly linear elasticity ... " << std::flush;
68
69 // Hole die Dichte \rho (density) und die Paramter \nu (Poisson-ratio) und \mu (zweite Lamé-Konstante)
70 double density = this->parameterList_->sublist("Parameter").get("Density",1000.);
71
72 double poissonRatio = this->parameterList_->sublist("Parameter").get("Poisson Ratio",0.4);
73 double mu = this->parameterList_->sublist("Parameter").get("Mu",2.0e+6);
74
75 // Berechne daraus nun E (Youngsches Modul) und die erste Lamé-Konstanten \lambda
76 double youngModulus = mu*2.*(1 + poissonRatio);
77 double lambda = (poissonRatio*youngModulus)/((1 + poissonRatio)*(1 - 2*poissonRatio));
78
79 // Initialisiere die Steifigkeitsmatrix. Das letzte Argument gibt die (ungefaehre) Anzahl an Eintraege pro Zeile an
80 MatrixPtr_Type K = Teuchos::rcp(new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getDimension() * this->getDomain(0)->getApproxEntriesPerRow() ) );
81 // MatrixPtr_Type K = Teuchos::rcp(new Matrix_Type( this->domainPtr_vec_.at(0)->getMapVecFieldUnique(), 10 ) );
82
83 // Assembliere die Steifigkeitsmatrix. Die 2 gibt degree an, d.h. die Ordnung der Quadraturformel, die benutzt werden soll.
84 this->feFactory_->assemblyLinElasXDim( this->dim_, this->getDomain(0)->getFEType(), K, lambda, mu );
85
86 // Setup fuer die linke Seite des zu loesdenen GLS. Beachte, dass system_ via system_() im Standardkonstruktor (von der Klasse Problem)
87 // initialisiert worden ist, hier wird dieser also erst richtig initialisiert.
88 // Ein Objekt der Klasse Bmat ist eine Blockmatrix; also ist system_ eine Blockmatrix (Objekt von BMat)
89
90 // Fuege die Steifikeitsmatrix als Blockeintrag an der Stelle (1,1) (in C dann (0,0)) in die Blockmatrix hinein.
91 this->system_->addBlock( K, 0, 0 );
92
93 this->assembleSourceTerm( 0. );
94 //this->sourceTerm_->scale(density);
95 this->addToRhs( this->sourceTerm_ );
96
97 if (this->verbose_)
98 std::cout << "done -- " << std::endl;
99}
100
101//template<class SC,class LO,class GO,class NO>
102//void LinElas<SC,LO,GO,NO>::assembleSourceTerm(double time)
103//{
104//
105// double force = this->parameterList_->sublist("Parameter").get("Volume force",0.);
106//
107// MultiVectorPtr_Type sourceTerm = Teuchos::rcp(new MultiVector_Type( this->domainPtr_vec_.at(0)->getMapVecFieldRepeated() ) );
108// vec_dbl_Type funcParameter(3,0.);
109//
110// funcParameter[0] = 0.; // degree of function
111// funcParameter[1] = time;
112// funcParameter[2] = force;
113//
114// this->feFactory_->assemblyRHS(this->dim_, this->domain_FEType_vec_.at(0), sourceTerm, "Vector", this->rhsFuncVec_[0], funcParameter);
115//
116// this->sourceTerm_->getBlockNonConst(0)->putScalar(0.);
117//
118// this->sourceTerm_->getBlockNonConst(0)->exportFromVector( sourceTerm, false, "Add" );
119//
120// double density = this->parameterList_->sublist("Parameter").get("Density",1.);
121//
122// this->sourceTerm_->getBlockNonConst(0)->scale(density);
123// // Noch mit \rho = density skalieren
124// // feRhs->scale(density);
125// // BlockMultiVectorPtr_Type blockMVSourceT = Teuchos::rcp( new BlockMultiVector_Type( 1 ) );
126// // blockMVSourceT->addBlock( feRhs, 0 );
127//
128// // Setze die rechte Seite auf das Attribut des SourceTerms aus dem Problem
129// // Dieses Attribut wird in der Zeitintegration benoetigt, vgl. DAESolverInTime
130// // this->sourceTerm_ = blockMVSourceT;
131//
140//}
141
142}
143#endif
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement.cpp:5