1#ifndef GEOMETRY_def_hpp
2#define GEOMETRY_def_hpp
3#include "Geometry_decl.hpp"
6void ZeroFErhsFunc2D(
double* x,
double* result,
double* parameters)
15void ZeroFErhsFunc3D(
double* x,
double* result,
double* parameters)
24double HeuristicScaling(
double* x,
double* parameter)
26 if(x[0] < parameter[0])
37template<
class SC,
class LO,
class GO,
class NO>
38Geometry<SC,LO,GO,NO>::Geometry(
const DomainConstPtr_Type &domain, std::string FEType, ParameterListPtr_Type parameterList):
39Problem_Type(parameterList,domain->getComm())
41 this->addVariable( domain , FEType ,
"d_f" , domain->getDimension());
42 this->dim_ = this->getDomain(0)->getDimension();
46template<
class SC,
class LO,
class GO,
class NO>
47Geometry<SC,LO,GO,NO>::~Geometry()
52template<
class SC,
class LO,
class GO,
class NO>
53void Geometry<SC,LO,GO,NO>::info(){
57template<
class SC,
class LO,
class GO,
class NO>
58void Geometry<SC,LO,GO,NO>::assemble( std::string type )
const
61 MatrixPtr_Type H = Teuchos::rcp(
new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getDimension() * this->getDomain(0)->getApproxEntriesPerRow() ) );
63 if (this->parameterList_->sublist(
"Parameter").get(
"Model",
"Laplace")==
"Laplace"){
65 std::cout <<
"-- Assembly Geometry (scaled Laplace)... " << std::flush;
66 double distanceLaplace = this->parameterList_->sublist(
"Parameter").get(
"Distance Laplace",0.03);
67 double coefficientLaplace = this->parameterList_->sublist(
"Parameter").get(
"Coefficient Laplace",1000.);
69 std::cout <<
"\n Distance Laplace = " << distanceLaplace <<
" coefficient Laplace = " << coefficientLaplace <<
" ... " << std::flush;
71 vec_dbl_Type parameter(2);
72 parameter[0] = distanceLaplace;
73 parameter[1] = coefficientLaplace;
75 this->feFactory_->assemblyLaplaceXDim(this->dim_, this->getDomain(0)->getFEType(), H, HeuristicScaling, &(parameter.at(0)) );
77 else if (this->parameterList_->sublist(
"Parameter").get(
"Model",
"Laplace")==
"Elasticity"){
79 std::cout <<
"-- Assembly Geometry (linear elasticity)... " << std::flush;
80 double poissonRatio = this->parameterList_->sublist(
"Parameter").get(
"Poisson Ratio",0.3);
81 double mu = this->parameterList_->sublist(
"Parameter").get(
"Mu",2.0e+6);
82 double youngModulus = mu*2.*(1 + poissonRatio);
83 double lambda = (poissonRatio*youngModulus)/((1 + poissonRatio)*(1 - 2*poissonRatio));
85 std::cout <<
"\n Poisson ratio = " << poissonRatio <<
" mu = "<<mu <<
" lambda = "<< lambda <<
" E = " << youngModulus <<
"... " << std::flush;
88 this->feFactory_->assemblyLinElasXDim( this->dim_, this->getDomain(0)->getFEType(), H, lambda, mu );
91 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Unknown model for geometry.");
93 this->system_->addBlock( H, 0, 0 );
95 this->assembleSourceTerm( 0. );
97 this->addToRhs( this->sourceTerm_ );
100 std::cout <<
"done -- " << std::endl;
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement.cpp:5