00001 #ifndef DUNE_AMGTRANSFER_HH
00002 #define DUNE_AMGTRANSFER_HH
00003
00004 #include<dune/istl/bvector.hh>
00005 #include<dune/istl/matrixredistribute.hh>
00006 #include<dune/istl/paamg/pinfo.hh>
00007 #include<dune/istl/owneroverlapcopy.hh>
00008 #include<dune/istl/paamg/aggregates.hh>
00009 #include<dune/common/exceptions.hh>
00010
00011 namespace Dune
00012 {
00013 namespace Amg
00014 {
00015
00026 template<class V1, class V2, class T>
00027 class Transfer
00028 {
00029
00030 public:
00031 typedef V1 Vertex;
00032 typedef V2 Vector;
00033
00034 template<typename T1, typename R>
00035 static void prolongate(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
00036 Vector& fineRedist,T1 damp, R& redistributor=R());
00037
00038 template<typename T1, typename R>
00039 static void prolongate(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
00040 T1 damp);
00041
00042 static void restrict(const AggregatesMap<Vertex>& aggregates, Vector& coarse, const Vector& fine,
00043 T& comm);
00044 };
00045
00046 template<class V,class V1>
00047 class Transfer<V,V1, SequentialInformation>
00048 {
00049 public:
00050 typedef V Vertex;
00051 typedef V1 Vector;
00052 typedef RedistributeInformation<SequentialInformation> Redist;
00053 template<typename T1>
00054 static void prolongate(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
00055 Vector& fineRedist, T1 damp,
00056 const SequentialInformation& comm=SequentialInformation(),
00057 const Redist& redist=Redist());
00058 template<typename T1>
00059 static void prolongate(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
00060 T1 damp,
00061 const SequentialInformation& comm=SequentialInformation());
00062
00063
00064 static void restrict(const AggregatesMap<Vertex>& aggregates, Vector& coarse, const Vector& fine,
00065 const SequentialInformation& comm);
00066 };
00067
00068 #if HAVE_MPI
00069
00070 template<class V,class V1, class T1, class T2>
00071 class Transfer<V,V1,OwnerOverlapCopyCommunication<T1,T2> >
00072 {
00073 public:
00074 typedef V Vertex;
00075 typedef V1 Vector;
00076 typedef RedistributeInformation<OwnerOverlapCopyCommunication<T1,T2> > Redist;
00077 template<typename T3>
00078 static void prolongate(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
00079 Vector& fineRedist, T3 damp, OwnerOverlapCopyCommunication<T1,T2>& comm,
00080 const Redist& redist);
00081 template<typename T3>
00082 static void prolongate(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
00083 T3 damp, OwnerOverlapCopyCommunication<T1,T2>& comm);
00084
00085 static void restrict(const AggregatesMap<Vertex>& aggregates, Vector& coarse, const Vector& fine,
00086 OwnerOverlapCopyCommunication<T1,T2>& comm);
00087 };
00088
00089 #endif
00090
00091 template<class V, class V1>
00092 template<typename T>
00093 inline void
00094 Transfer<V,V1,SequentialInformation>::prolongate(const AggregatesMap<Vertex>& aggregates,
00095 Vector& coarse, Vector& fine, Vector& fineRedist,
00096 T damp,
00097 const SequentialInformation& comm,
00098 const Redist& redist)
00099 {
00100 prolongate(aggregates, coarse, fine, damp);
00101 }
00102 template<class V, class V1>
00103 template<typename T>
00104 inline void
00105 Transfer<V,V1,SequentialInformation>::prolongate(const AggregatesMap<Vertex>& aggregates,
00106 Vector& coarse, Vector& fine,
00107 T damp,
00108 const SequentialInformation& comm)
00109 {
00110 typedef typename Vector::iterator Iterator;
00111
00112 Iterator end = coarse.end();
00113 Iterator begin= coarse.begin();
00114 for(;begin!=end;++begin)
00115 *begin*=damp;
00116 end=fine.end();
00117 begin=fine.begin();
00118
00119 for(Iterator block=begin; block != end; ++block){
00120 std::ptrdiff_t index=block-begin;
00121 const Vertex& vertex = aggregates[index];
00122 if(vertex != AggregatesMap<Vertex>::ISOLATED)
00123 *block += coarse[aggregates[index]];
00124 }
00125 }
00126
00127 template<class V, class V1>
00128 inline void
00129 Transfer<V,V1,SequentialInformation>::restrict(const AggregatesMap<Vertex>& aggregates,
00130 Vector& coarse,
00131 const Vector& fine,
00132 const SequentialInformation& comm)
00133 {
00134
00135 coarse=0;
00136
00137 typedef typename Vector::const_iterator Iterator;
00138 Iterator end = fine.end();
00139 Iterator begin=fine.begin();
00140
00141 for(Iterator block=begin; block != end; ++block){
00142 const Vertex& vertex = aggregates[block-begin];
00143 if(vertex != AggregatesMap<Vertex>::ISOLATED)
00144 coarse[vertex] += *block;
00145 }
00146 }
00147
00148 #if HAVE_MPI
00149 template<class V, class V1, class T1, class T2>
00150 template<typename T3>
00151 inline void Transfer<V,V1,OwnerOverlapCopyCommunication<T1,T2> >::prolongate(const AggregatesMap<Vertex>& aggregates,
00152 Vector& coarse, Vector& fine,
00153 Vector& fineRedist, T3 damp,
00154 OwnerOverlapCopyCommunication<T1,T2>& comm,
00155 const Redist& redist)
00156 {
00157 if(fineRedist.size()>0)
00158
00159 Transfer<V,V1,SequentialInformation>::prolongate(aggregates, coarse, fineRedist, damp);
00160
00161
00162 redist.redistributeBackward(fine, fineRedist);
00163 comm.copyOwnerToAll(fine,fine);
00164 }
00165
00166 template<class V, class V1, class T1, class T2>
00167 template<typename T3>
00168 inline void Transfer<V,V1,OwnerOverlapCopyCommunication<T1,T2> >::prolongate(const AggregatesMap<Vertex>& aggregates,
00169 Vector& coarse, Vector& fine,
00170 T3 damp,
00171 OwnerOverlapCopyCommunication<T1,T2>& comm)
00172 {
00173 Transfer<V,V1,SequentialInformation>::prolongate(aggregates, coarse, fine, damp);
00174 }
00175 template<class V, class V1, class T1, class T2>
00176 inline void Transfer<V,V1,OwnerOverlapCopyCommunication<T1,T2> >::restrict(const AggregatesMap<Vertex>& aggregates,
00177 Vector& coarse, const Vector& fine,
00178 OwnerOverlapCopyCommunication<T1,T2>& comm)
00179 {
00180 Transfer<V,V1,SequentialInformation>::restrict(aggregates, coarse, fine, SequentialInformation());
00181
00182
00183 comm.project(coarse);
00184 }
00185 #endif
00186
00187 }
00188 }
00189 #endif