Loading [MathJax]/extensions/tex2jax.js

dune-composites (unstable)

geneo.hh
1// -*- tab-width: 4; indent-tabs-mode: nil -*-
2// vi: set et ts=4 sw=4 sts=4:
3#ifndef GENEO_SOLVER_HH
4#define GENEO_SOLVER_HH
5
6#include <dune/pdelab/backend/istl/geneo/geneo.hh>
7#include "zembasis.hh"
8
9namespace Dune{
10 namespace Geneo{
11
12 template<class V, class GFS, class GFS_EXT, class GO, class GO_EXT, class LOP, class CON, class MBE>
13 class CG_GenEO :
14 public PDELab::OVLPScalarProductImplementation<GFS>,
15 public PDELab::LinearResultStorage
16 {
17
18 public:
19
20 typedef double RF;
21
22 CG_GenEO(V& x0_, const GFS& gfs_, const GFS_EXT gfs_ext_, LOP& lop_, CON& constraints_, MBE& mbe_, Composites::SolverInfo & si_, Dune::MPIHelper& helper_, int cells_) :
23 PDELab::OVLPScalarProductImplementation<GFS>(gfs_),
24 x0(x0_),
25 gfs(gfs_),
26 gfs_ext(gfs_ext_),
27 lop(lop_),
28 constraints(constraints_),
29 mbe(mbe_),
30 solver_info(si_),
31 helper(helper_),
32 cells(cells_) { }
33
34 inline void apply()
35 {
36#if HAVE_SUITESPARSE
37 typedef typename GFS::template ConstraintsContainer<RF>::Type CC;
38 auto cg = CC();
39 typedef typename GFS_EXT::template ConstraintsContainer<RF>::Type CC_EXTERIOR;
40 auto cg_ext = CC_EXTERIOR();
41 auto cc_bnd_neu_int_dir = CC();
42
43 // assemble constraints
44 Dune::PDELab::constraints(constraints,gfs,cg);
45 Dune::PDELab::constraints(constraints,gfs_ext,cg_ext);
46
47 Dune::PDELab::NoDirichletConstraintsParameters pnbc;
48 Dune::PDELab::constraints(pnbc,gfs,cc_bnd_neu_int_dir);
49
50 std::vector<double> times;
51 Dune::Timer watch;
52 watch.reset();
53
54 // also create a wrapper using the same data based on the GFS without processor boundary constraints
55 typedef Dune::PDELab::Backend::Vector<GFS_EXT, RF> V_EXTERIOR;
56 V_EXTERIOR x0_ext(gfs_ext, Dune::PDELab::Backend::unattached_container());
57 x0_ext.attach(x0.storage());
58
59 //Set up Grid operator
60 typedef typename GO::Jacobian J;
61 typedef typename GO_EXT::Jacobian J_ext;
62
63 //Matrix with "correct" boundary conditions
64 GO go(gfs,cg,gfs,cg,lop,mbe);
65 J AF(go);
66 AF = 0.0;
67 go.jacobian(x0,AF);
68
69 //Matrix with pure neumann boundary conditions
70 GO_EXT go_neu(gfs_ext,cg_ext,gfs_ext,cg_ext,lop,mbe);
71 J_ext AF_neu(go_neu);
72 AF_neu = 0.0;
73 go_neu.jacobian(x0_ext,AF_neu);
74
75 //Create local operator on overlap
76 typedef Dune::PDELab::LocalOperatorOvlpRegion<LOP, GFS_EXT> LOP_ovlp;
77 LOP_ovlp lop_ovlp(lop,gfs_ext);
78 typedef Dune::PDELab::GridOperator<GFS_EXT,GFS_EXT,LOP_ovlp,MBE,RF,RF,RF,CC,CC> GO_ovlp;
79 GO_ovlp go_ovlp(gfs_ext,cg_ext,gfs_ext,cg_ext,lop_ovlp,mbe);
80 J_ext AF_ovlp(go_ovlp);
81 AF_ovlp = 0.0;
82 go_ovlp.jacobian(x0_ext,AF_ovlp);
83
84 //Set up solution vector and some necessary operators
85 typedef Dune::PDELab::OverlappingOperator<CC,J,V,V> POP;
86 POP popf(cg,AF);
87 typedef Dune::PDELab::ISTL::ParallelHelper<GFS> PIH;
88 PIH pihf(gfs);
89 typedef Dune::PDELab::OverlappingScalarProduct<GFS,V> OSP;
90 OSP ospf(gfs,pihf);
91
92 // auto gv = gfs.gridView();
93
94 double eigenvalue_threshold = solver_info.eigenvalue_threshold*(double)solver_info.overlap/(cells+solver_info.overlap); //eigenvalue threshhold
95 if(helper.rank() == 0){
96 if(eigenvalue_threshold > 0) std::cout << "Eigenvalue threshhold: " << eigenvalue_threshold << std::endl;
97 else std::cout << "Using " << solver_info.nev << " EV in each subdomain" << std::endl;
98 }
99 typedef Dune::PDELab::LocalFunctionSpace<GFS, Dune::PDELab::AnySpaceTag> LFSU;
100 typedef typename LFSU::template Child<0>::Type LFS;
101 LFSU lfsu(gfs);
102 LFS lfs = lfsu.template child<0>();
103
104 std::shared_ptr<V> part_unity;
105 part_unity = std::make_shared<V>(standardPartitionOfUnity<V>(gfs, cc_bnd_neu_int_dir));
106
107 V dummy(gfs, 1);
108 Dune::PDELab::set_constrained_dofs(cg_ext,0.0,dummy); // Zero on subdomain boundary
109
110 int nev = solver_info.nev;
111
112 std::shared_ptr<Dune::PDELab::SubdomainBasis<V> > subdomain_basis;
113 if (nev > 0)
114 subdomain_basis = std::make_shared<Dune::PDELab::GenEOBasis<GFS,J_ext,V,3> >(gfs, AF_neu, AF_ovlp, eigenvalue_threshold, *part_unity, solver_info.nev, solver_info.nev_arpack, 0.001, false, solver_info.verb);//,subdomain_basis2);
115 else if (nev == 0){
116 subdomain_basis = std::make_shared<Dune::PDELab::ZEMBasis<GFS,LFS,V,3,3> >(gfs, lfs, *part_unity);
117 }
118 else{
119 subdomain_basis = std::make_shared<Dune::PDELab::SubdomainBasis<V> >(*part_unity);
120 }
121
122 auto partunityspace = std::make_shared<Dune::PDELab::SubdomainProjectedCoarseSpace<GFS,J_ext,V,PIH> >(gfs, AF_neu, subdomain_basis, pihf);
123 typedef Dune::PDELab::ISTL::TwoLevelOverlappingAdditiveSchwarz<GFS,J,V,V> PREC_PCG;
124 std::shared_ptr<PREC_PCG> prec;
125 prec = std::make_shared<PREC_PCG>(gfs, AF, partunityspace,solver_info.coarseSpaceActive);
126
127 MPI_Barrier(gfs.gridView().comm());
128 if(helper.rank()==0) std::cout << "Geneo setup time " << watch.elapsed() << std::endl;
129 times.push_back(watch.elapsed());
130 watch.reset();
131
132 // set up and assemble right hand side w.r.t. l(v)-a(u_g,v)
133 V d(gfs,0.0);
134 go.residual(x0,d);
135 V v(gfs,0.0);
136
137 //Solve using CG
138 std::shared_ptr<CGSolver<V>> solver;
139 solver = std::make_shared<CGSolver<V> >(popf,ospf,*prec,solver_info.KrylovTol,solver_info.MaxIt,(helper.rank()==0) ? solver_info.verb : 0,true);
140 Dune::InverseOperatorResult result;
141 solver->apply(v,d,result);
142
143 times.push_back(watch.elapsed());
144 watch.reset();
145 x0 -= v;
146
147 solver_info.setTimes(times);
148 solver_info.recordResult(result);
149
150#else
151 std::cout << "Solver not available . . . Please install UMFPACK as part of SuiteSparse" << std::endl;
152 return;
153#endif
154 }
155
156 //Templated apply for Newton solver
157 template<class M, class V2, class W>
158 void apply(M& A, V2& z, W& r, typename Dune::template FieldTraits<typename V2::ElementType >::real_type reduction)
159 {
160 this->apply(z,r,reduction);
161 }
162
163 private:
164 V & x0;
165 const GFS& gfs;
166 const GFS_EXT& gfs_ext;
167 LOP& lop;
168 const CON& constraints;
169 const MBE& mbe;
170 Composites::SolverInfo& solver_info;
171 Dune::MPIHelper& helper;
172 int cells;
173 };
174 }
175}
176
177#endif
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Apr 5, 23:02, 2025)