00001 #ifdef CH_LANG_CC
00002
00003
00004
00005
00006
00007
00008
00009 #endif
00010
00011
00012
00013
00014
00015 #ifndef _LEVELDATAOPS_H_
00016 #define _LEVELDATAOPS_H_
00017
00018 #ifdef CH_MPI
00019 #include <string>
00020 #include <sstream>
00021 #endif
00022 #include "LevelData.H"
00023 #include "RefCountedPtr.H"
00024 #include "SPMD.H"
00025 #include "Copier.H"
00026 #include "NamespaceHeader.H"
00027
00028
00029
00030 template <class T>
00031 class LevelDataOps
00032 {
00033 public:
00034 LevelDataOps()
00035 :m_levelFactory( new DefaultDataFactory<T>() )
00036 {
00037 }
00038
00039 LevelDataOps(RefCountedPtr<DataFactory<T> > a_factoryPtr)
00040 :m_levelFactory(a_factoryPtr)
00041 {
00042 }
00043
00044 virtual ~LevelDataOps()
00045 {
00046 }
00047
00048 virtual void define(const RefCountedPtr<DataFactory<T> >& a_factoryPtr)
00049 {
00050 m_levelFactory = a_factoryPtr;
00051 }
00052
00053 virtual void define(DataFactory<T>* a_rawPointer)
00054 {
00055 m_levelFactory = RefCountedPtr<DataFactory<T> >(a_rawPointer);
00056 }
00057
00058 virtual void create( LevelData<T>& a_lhs, const LevelData<T>& a_rhs);
00059
00060 virtual void assign( LevelData<T>& a_lhs, const LevelData<T>& a_rhs) ;
00061
00062 virtual Real dotProduct(const LevelData<T>& a_1, const LevelData<T>& a_2) ;
00063
00064 virtual void mDotProduct(const LevelData<T>& a_1, const int a_sz, const LevelData<T> a_2arr[], Real a_mdots[]);
00065
00066 virtual void incr( LevelData<T>& a_lhs, const LevelData<T>& a_x, Real a_scale) ;
00067
00068 virtual void mult( LevelData<T>& a_lhs, const LevelData<T>& a_x);
00069
00070 virtual void axby( LevelData<T>& a_lhs, const LevelData<T>& a_x,
00071 const LevelData<T>& a_y, Real a_a, Real a_b) ;
00072
00073 virtual void scale(LevelData<T>& a_lhs, const Real& a_scale) ;
00074
00075 virtual void plus(LevelData<T>& a_lhs, const Real& a_inc) ;
00076
00077 virtual void setToZero(LevelData<T>& a_lhs);
00078
00079 virtual void setVal(LevelData<T>& a_lhs, const Real& a_val);
00080
00081 virtual void copyToZero(LevelData<T>& a_lhs, const Copier& a_copier);
00082
00083
00084 protected:
00085 RefCountedPtr<DataFactory<T> > m_levelFactory;
00086 };
00087
00088
00089
00090
00091
00092
00093 template <class T>
00094 void LevelDataOps<T>:: create(LevelData<T>& a_lhs, const LevelData<T>& a_rhs)
00095 {
00096
00097 a_lhs.define(a_rhs.disjointBoxLayout(), a_rhs.nComp(),
00098 a_rhs.ghostVect(), *m_levelFactory);
00099 }
00100
00101 template <class T>
00102 void LevelDataOps<T>:: assign(LevelData<T>& a_lhs, const LevelData<T>& a_rhs)
00103 {
00104 Interval interv(0, a_rhs.nComp()-1);
00105 a_rhs.copyTo(interv, a_lhs, interv);
00106 }
00107
00108 template <class T>
00109 Real LevelDataOps<T>::dotProduct(const LevelData<T>& a_1, const LevelData<T>& a_2)
00110 {
00111 const DisjointBoxLayout& dbl = a_1.disjointBoxLayout();
00112 Real val = 0.0;
00113 DataIterator dit=dbl.dataIterator(); int ompsize=dit.size();
00114 #pragma omp parallel for reduction (+:val)
00115 for(int i=0; i<ompsize; i++)
00116 {
00117 const DataIndex& d = dit[i];
00118 val += a_1[d].dotProduct(a_2[d], dbl.get(d));
00119 }
00120
00121
00122 #ifdef CH_MPI
00123 Real recv;
00124 int result = MPI_Allreduce(&val, &recv, 1, MPI_CH_REAL,
00125 MPI_SUM, Chombo_MPI::comm);
00126 if ( result != 0 )
00127 {
00128 std::ostringstream msg;
00129 msg << "LevelDataOps::dotProduct() called MPI_Allreduce() which returned error code " << result ;
00130 MayDay::Warning( msg.str().c_str() );
00131 }
00132 val = recv;
00133 #endif
00134 return val;
00135
00136 }
00137
00138
00139 template <class T>
00140 void LevelDataOps<T>::mDotProduct(const LevelData<T>& a_1, const int a_sz, const LevelData<T> a_2arr[], Real a_mdots[])
00141 {
00142 const DisjointBoxLayout& dbl = a_1.disjointBoxLayout();
00143
00144 for (int ii=0;ii<a_sz;ii++)
00145 {
00146 Real val = 0.0;
00147 const LevelData<T> &a_2 = a_2arr[ii];
00148 DataIterator dit=dbl.dataIterator(); int ompsize=dit.size();
00149 #pragma omp parallel for reduction (+:val)
00150 for(int i=0; i<ompsize; i++)
00151 {
00152 const DataIndex& d = dit[i];
00153 val += a_1[d].dotProduct(a_2[d], dbl.get(d));
00154 }
00155 a_mdots[ii] = val;
00156 }
00157
00158 #ifdef CH_MPI
00159 Real *recv = new Real[a_sz];
00160
00161 int result = MPI_Allreduce(a_mdots, recv, a_sz, MPI_CH_REAL,
00162 MPI_SUM, Chombo_MPI::comm);
00163 if ( result != 0 )
00164 {
00165 std::ostringstream msg;
00166 msg << "LevelDataOps::mDotProduct() called MPI_Allreduce() which returned error code " << result ;
00167 MayDay::Warning( msg.str().c_str() );
00168 }
00169 for (int ii=0;ii<a_sz;ii++)
00170 {
00171 a_mdots[ii] = recv[ii];
00172 }
00173
00174 delete [] recv;
00175 #endif
00176 }
00177
00178 template <class T>
00179 void LevelDataOps<T>:: incr( LevelData<T>& a_lhs, const LevelData<T>& a_rhs, Real a_scale)
00180 {
00181
00182 int numcomp = a_lhs.nComp();
00183 int startcomp = 0;
00184 DataIterator dit=a_lhs.dataIterator(); int count=dit.size();
00185 #pragma omp parallel for
00186 for(int i=0; i<count; i++)
00187 {
00188 const DataIndex& d=dit[i];
00189 Box subbox = a_lhs.disjointBoxLayout()[d];
00190 a_lhs[d].plus(a_rhs[d], subbox, subbox, a_scale, startcomp, startcomp, numcomp);
00191 }
00192 }
00193
00194 template <class T>
00195 void LevelDataOps<T>:: mult( LevelData<T>& a_lhs, const LevelData<T>& a_rhs)
00196 {
00197
00198 DataIterator dit=a_lhs.dataIterator(); int count=dit.size();
00199 #pragma omp parallel for
00200 for(int i=0; i<count; i++)
00201 {
00202 const DataIndex& d=dit[i];
00203 a_lhs[d] *= a_rhs[d];
00204 }
00205 }
00206
00207 template <class T>
00208 void LevelDataOps<T>:: setToZero( LevelData<T>& a_lhs)
00209 {
00210 DataIterator dit=a_lhs.dataIterator(); int count=dit.size();
00211 #pragma omp parallel for
00212 for(int i=0; i<count; i++)
00213 {
00214 const DataIndex& d=dit[i];
00215 a_lhs[d].setVal(0.0);
00216 }
00217 }
00218 template <class T>
00219 void LevelDataOps<T>:: setVal( LevelData<T>& a_lhs, const Real& a_val)
00220 {
00221 DataIterator dit=a_lhs.dataIterator(); int count=dit.size();
00222 #pragma omp parallel for
00223 for(int i=0; i<count; i++)
00224 {
00225 const DataIndex& d=dit[i];
00226 a_lhs[d].setVal(a_val);
00227 }
00228 }
00229
00230 template <class T>
00231 void LevelDataOps<T>:: copyToZero( LevelData<T>& a_lhs, const Copier& a_copier)
00232 {
00233 int nComp = a_lhs.nComp();
00234 CopyIterator cit(a_copier, CopyIterator::LOCAL); int ccount=cit.size();
00235 CopyIterator tit(a_copier, CopyIterator::TO); int tcount = tit.size();
00236
00237 #pragma omp parallel
00238 {
00239 #pragma omp for nowait
00240 for (int i=0; i<ccount; i++)
00241 {
00242 const MotionItem& item = cit[i];
00243 a_lhs[item.toIndex].setVal(0, item.toRegion, 0, nComp);
00244 }
00245 #pragma omp for nowait
00246 for (int i=0; i<tcount; i++)
00247 {
00248 const MotionItem& item = tit[i];
00249 a_lhs[item.toIndex].setVal(0, item.toRegion, 0, nComp);
00250 }
00251 }
00252 }
00253
00254 template <class T>
00255 void LevelDataOps<T>:: axby( LevelData<T>& a_lhs, const LevelData<T>& a_x,
00256 const LevelData<T>& a_y, Real a, Real b)
00257 {
00258 DataIterator dit=a_lhs.dataIterator(); int count=dit.size();
00259 #pragma omp parallel for
00260 for(int i=0; i<count; i++)
00261 {
00262 const DataIndex& d=dit[i];
00263 T& data = a_lhs[d];
00264 data.copy(a_x[d]);
00265 data.mult(a);
00266 data.plus(a_y[d], b);
00267 }
00268 }
00269
00270 template <class T>
00271 void LevelDataOps<T>:: scale(LevelData<T>& a_lhs, const Real& a_scale)
00272 {
00273 DataIterator dit=a_lhs.dataIterator(); int count=dit.size();
00274 #pragma omp parallel for
00275 for(int i=0; i<count; i++)
00276 {
00277 const DataIndex& d=dit[i];
00278 T& data = a_lhs[d];
00279 data.mult(a_scale);
00280 }
00281 }
00282
00283 template <class T>
00284 void LevelDataOps<T>:: plus(LevelData<T>& a_lhs, const Real& a_inc)
00285 {
00286 DataIterator dit=a_lhs.dataIterator(); int count=dit.size();
00287 #pragma omp parallel for
00288 for(int i=0; i<count; i++)
00289 {
00290 const DataIndex& d=dit[i];
00291 T& data = a_lhs[d];
00292 data += a_inc;
00293 }
00294 }
00295
00296
00297
00298 #include "NamespaceFooter.H"
00299 #endif