Loading [MathJax]/extensions/tex2jax.js

dune-composites (unstable)

plot_properties.hh
1#ifndef PLOT_STRESSES_HH
2#define PLOT_STRESSES_HH
3
4#include "composite_properties_lop.hh"
5
6namespace Dune{
7 namespace Composites{
8
9 template<class MODEL, typename V, typename GV1, typename GFS, typename MBE>
10 void plotProperties(MODEL& model, const GV1& gv1, const GFS& gfs, MBE& mbe){
11 using Dune::PDELab::Backend::native;
12
13 typedef Dune::PDELab::ISTL::VectorBackend<> Scalar_VectorBackend;
14 typedef double RF;
15 std::string vtk_output = config.get<std::string>("ModelName","TestModel");
16 Dune::Timer timer;
17
18 using GV = Dune::PDELab::AllEntitySet<GV1>;
19 auto gv = GV(gv1);
20
21 typedef typename GV::Grid::ctype e_ctype;
22 typedef Dune::PDELab::QkDGLocalFiniteElementMap<e_ctype,double,0,3> FEM_properties;
23 FEM_properties fem_properties;
24
25 typedef Dune::PDELab::GridFunctionSpace<GV, FEM_properties, Dune::PDELab::NoConstraints, Scalar_VectorBackend> SCALAR_P1GFS;
26 SCALAR_P1GFS materialType(gv,fem_properties); materialType.name("materialType");
27 SCALAR_P1GFS orientation(gv,fem_properties); orientation.name("orientation");
28
29 typedef Dune::PDELab::ISTL::VectorBackend<Dune::PDELab::ISTL::Blocking::fixed> VectorBackend;
30
31 typedef Dune::PDELab::PowerGridFunctionSpace <SCALAR_P1GFS, 2, VectorBackend,Dune::PDELab::EntityBlockedOrderingTag> P1GFS;
32
33 P1GFS p1gfs(materialType,orientation);
34
35 typedef Dune::PDELab::EmptyTransformation NoTrafo;
36
37 Dune::PDELab::composite_properties<GV,MODEL> lopProperties(gv,model);
38
39 typedef Dune::PDELab::GridOperator<P1GFS,P1GFS,Dune::PDELab::composite_properties<GV,MODEL>,MBE,RF,RF,RF,NoTrafo,NoTrafo> GO;
40 GO go(p1gfs,p1gfs,lopProperties,mbe);
41
42 typedef typename GO::Traits::Range V1;
43 V1 properties(p1gfs,0.0), dummy(p1gfs,0.0);
44 go.residual(dummy,properties);
45
46 Dune::SubsamplingVTKWriter<GV1> vtkwriter(gv1,Dune::refinementLevels(0));
47 Dune::PDELab::addSolutionToVTKWriter(vtkwriter,p1gfs,properties);
48 vtkwriter.write(model.vtk_properties_output,Dune::VTK::appendedraw);
49 }
50
51
52 template<class MODEL, typename V, typename GV1, typename GFS, typename MBE>
53 void plotPropertiesUG(MODEL& model, const GV1& gv1, const GFS& gfs, MBE& mbe){
54 using Dune::PDELab::Backend::native;
55
56 typedef Dune::PDELab::ISTL::VectorBackend<> Scalar_VectorBackend;
57 typedef double RF;
58 std::string vtk_output = config.get<std::string>("ModelName","TestModel");
59 Dune::Timer timer;
60
61 using GV = Dune::PDELab::AllEntitySet<GV1>;
62 auto gv = GV(gv1);
63
64 typedef typename GV::Grid::ctype e_ctype;
65 typedef Dune::PDELab::PkLocalFiniteElementMap<GV,e_ctype,double,1> FEM_properties;
66 FEM_properties fem_properties(gv);
67
68 typedef Dune::PDELab::GridFunctionSpace<GV, FEM_properties, Dune::PDELab::NoConstraints, Scalar_VectorBackend> SCALAR_P1GFS;
69 SCALAR_P1GFS materialType(gv,fem_properties); materialType.name("materialType");
70 SCALAR_P1GFS orientation(gv,fem_properties); orientation.name("orientation");
71
72 typedef Dune::PDELab::ISTL::VectorBackend<Dune::PDELab::ISTL::Blocking::fixed> VectorBackend;
73
74 typedef Dune::PDELab::PowerGridFunctionSpace <SCALAR_P1GFS, 2, VectorBackend,Dune::PDELab::EntityBlockedOrderingTag> P1GFS;
75
76 // Define operator to count elements containing a degree of freeedom
77 Dune::PDELab::countElem lopCnt;
78 SCALAR_P1GFS cntpgfs(gv,fem_properties);
79 typedef Dune::PDELab::EmptyTransformation NoTrafo;
80 typedef Dune::PDELab::GridOperator<GFS,SCALAR_P1GFS,Dune::PDELab::countElem,MBE,RF,RF,RF,NoTrafo,NoTrafo> GOCnt;
81 GOCnt goCnt(gfs,cntpgfs,lopCnt,mbe);
82
83 typedef typename GOCnt::Traits::Range Vcnt;
84 Vcnt ElemCnt(cntpgfs,0.0);
85
86 V x0(gfs,0.0);
87 goCnt.residual(x0,ElemCnt);
88
89 //Turn counter into discrete grid function
90 typedef Dune::PDELab::DiscreteGridFunction<SCALAR_P1GFS,Vcnt> DGF;
91 DGF dgfCnt(cntpgfs, ElemCnt);
92
93 // Set up plot properties
94 P1GFS p1gfs(materialType,orientation);
95
96 typedef Dune::PDELab::EmptyTransformation NoTrafo;
97
98 Dune::PDELab::composite_propertiesUG<GV,MODEL, DGF> lopProperties(gv,model,dgfCnt);
99
100 typedef Dune::PDELab::GridOperator<P1GFS,P1GFS,Dune::PDELab::composite_propertiesUG<GV,MODEL,DGF>,MBE,RF,RF,RF,NoTrafo,NoTrafo> GO;
101 GO go(p1gfs,p1gfs,lopProperties,mbe);
102
103 typedef typename GO::Traits::Range V1;
104 V1 properties(p1gfs,0.0), dummy(p1gfs,0.0);
105 go.residual(dummy,properties);
106
107 Dune::SubsamplingVTKWriter<GV1> vtkwriter(gv1,Dune::refinementLevels(0));
108 Dune::PDELab::addSolutionToVTKWriter(vtkwriter,p1gfs,properties);
109 vtkwriter.write(model.vtk_properties_output,Dune::VTK::appendedraw);
110 }
111
112 }
113}
114
115#endif
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Apr 3, 22:46, 2025)