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