• Main Page
  • Related Pages
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

Albany_Catalyst_TeuchosArrayRCPDataArray.hpp

Go to the documentation of this file.
00001 //*****************************************************************//
00002 //    Albany 2.0:  Copyright 2013 Kitware Inc.                     //
00003 //    This Software is released under the BSD license detailed     //
00004 //    in the file "license.txt" in the top-level Albany directory  //
00005 //*****************************************************************//
00006 
00007 #include "Albany_Catalyst_TeuchosArrayRCPDataArrayDecl.hpp"
00008 
00009 #include "Shards_BasicTopologies.hpp"
00010 
00011 #include "Teuchos_TestForException.hpp"
00012 
00013 #include "vtkIdList.h"
00014 #include "vtkIdTypeArray.h"
00015 #include "vtkObjectFactory.h"
00016 #include "vtkVariantCast.h"
00017 
00018 #include <cassert>
00019 #include <exception>
00020 
00021 namespace Albany {
00022 namespace Catalyst {
00023 
00024 // Can't use vtkStandardNewMacro on a templated class:
00025 template <typename Scalar>
00026 TeuchosArrayRCPDataArray<Scalar> * TeuchosArrayRCPDataArray<Scalar>::New()
00027 {
00028   VTK_STANDARD_NEW_BODY(TeuchosArrayRCPDataArray<Scalar>)
00029 }
00030 
00031 template <typename Scalar>
00032 void TeuchosArrayRCPDataArray<Scalar>::PrintSelf(ostream &os, vtkIndent indent)
00033 {
00034   this->TeuchosArrayRCPDataArray<Scalar>::Superclass::PrintSelf(os, indent);
00035 }
00036 
00037 template <typename Scalar>
00038 void TeuchosArrayRCPDataArray<Scalar>::
00039 SetArrayRCP(const Teuchos::ArrayRCP<Scalar> &array, int numComps)
00040 {
00041   this->Data = array;
00042   this->NumberOfComponents = numComps;
00043   delete [] this->TmpArray;
00044   this->TmpArray = new Scalar[numComps];
00045   this->Size = static_cast<vtkIdType>(array.size());
00046   this->MaxId = this->Size - 1;
00047 }
00048 
00049 template <typename Scalar>
00050 void TeuchosArrayRCPDataArray<Scalar>::Initialize()
00051 {
00052   this->Data = Teuchos::ArrayRCP<Scalar>();
00053   this->Size = 0;
00054   this->MaxId = -1;
00055 }
00056 
00057 template <typename Scalar>
00058 void TeuchosArrayRCPDataArray<Scalar>::GetTuples(vtkIdList *ptIds,
00059                                                  vtkAbstractArray *output)
00060 {
00061   vtkDataArray *da = vtkDataArray::FastDownCast(output);
00062   TEUCHOS_TEST_FOR_EXCEPTION(!da, std::runtime_error,
00063                              "GetTuples must be given a vtkDataArray arg.");
00064 
00065   da->Reset();
00066   da->SetNumberOfComponents(this->NumberOfComponents);
00067   da->Allocate(this->NumberOfComponents * ptIds->GetNumberOfIds());
00068 
00069   vtkIdType *begin = ptIds->GetPointer(0);
00070   vtkIdType *end = ptIds->GetPointer(ptIds->GetNumberOfIds());
00071 
00072   // Use type appropriate API if possible:
00073   if (vtkTypedDataArray<Scalar> *ta
00074       = vtkTypedDataArray<Scalar>::FastDownCast(da)) {
00075     while (begin != end) {
00076       ta->InsertNextTupleValue(
00077             &this->Data[*(begin++) * this->NumberOfComponents]);
00078     }
00079   }
00080   else { // otherwise, use the double interface:
00081     double *tuple = new double[this->NumberOfComponents];
00082     while (begin != end) {
00083       this->GetTuple(*(begin++), tuple);
00084       da->InsertNextTuple(tuple);
00085     }
00086     delete [] tuple;
00087   }
00088 }
00089 
00090 template <typename Scalar>
00091 void TeuchosArrayRCPDataArray<Scalar>::GetTuples(vtkIdType p1, vtkIdType p2,
00092                                                  vtkAbstractArray *output)
00093 {
00094   vtkDataArray *da = vtkDataArray::FastDownCast(output);
00095   TEUCHOS_TEST_FOR_EXCEPTION(!da, std::runtime_error,
00096                              "GetTuples must be given a vtkDataArray arg.");
00097 
00098   da->Reset();
00099   da->SetNumberOfComponents(this->NumberOfComponents);
00100   da->Allocate(this->NumberOfComponents * (p2 - p1 + 1));
00101 
00102   // Use type appropriate API if possible:
00103   if (vtkTypedDataArray<Scalar> *ta
00104       = vtkTypedDataArray<Scalar>::FastDownCast(da)) {
00105     while (p1 <= p2)
00106       ta->InsertNextTupleValue(&this->Data[p1++ * this->NumberOfComponents]);
00107   }
00108   else { // otherwise, use the double interface:
00109     double *tuple = new double[this->NumberOfComponents];
00110     while (p1 <= p2) {
00111       this->GetTuple(p1++, tuple);
00112       da->InsertNextTuple(tuple);
00113     }
00114     delete [] tuple;
00115   }
00116 }
00117 
00118 template <typename Scalar>
00119 void TeuchosArrayRCPDataArray<Scalar>::Squeeze()
00120 {
00121 }
00122 
00123 template <typename Scalar>
00124 vtkArrayIterator *TeuchosArrayRCPDataArray<Scalar>::NewIterator()
00125 {
00126   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error,
00127                              "Method not implemented.");
00128   return NULL;
00129 }
00130 
00131 template <typename Scalar>
00132 vtkIdType TeuchosArrayRCPDataArray<Scalar>::LookupValue(vtkVariant value)
00133 {
00134   bool valid = false;
00135   double tmp = vtkVariantCast<double>(value, &valid);
00136   if (valid)
00137     return this->Lookup(tmp, 0);
00138   return -1;
00139 }
00140 
00141 template <typename Scalar>
00142 void TeuchosArrayRCPDataArray<Scalar>::
00143 LookupValue(vtkVariant value, vtkIdList *ids)
00144 {
00145   bool valid = false;
00146   double tmp = vtkVariantCast<double>(value, &valid);
00147   ids->Reset();
00148   if (valid) {
00149     vtkIdType index = 0;
00150     while ((index = this->Lookup(tmp, index)) >= 0)
00151       ids->InsertNextId(index++);
00152   }
00153 }
00154 
00155 template <typename Scalar>
00156 vtkVariant TeuchosArrayRCPDataArray<Scalar>::GetVariantValue(vtkIdType idx)
00157 {
00158   return vtkVariant(this->Data[idx]);
00159 }
00160 
00161 template <typename Scalar>
00162 void TeuchosArrayRCPDataArray<Scalar>::ClearLookup()
00163 {
00164   // no-op, no fast lookup implemented.
00165 }
00166 
00167 template <typename Scalar>
00168 double *TeuchosArrayRCPDataArray<Scalar>::GetTuple(vtkIdType i)
00169 {
00170   int numComps = this->NumberOfComponents;
00171   double *out = this->TmpArray;
00172   typename ContainerType::const_iterator in =
00173       this->Data.begin() + (i * numComps);
00174   while (numComps-- > 0)
00175     *(out++) = static_cast<double>(*(in++));
00176 
00177   return this->TmpArray;
00178 }
00179 
00180 template <typename Scalar>
00181 void TeuchosArrayRCPDataArray<Scalar>::GetTuple(vtkIdType i, double *tuple)
00182 {
00183   int numComps = this->NumberOfComponents;
00184   typename ContainerType::const_iterator in
00185       = this->Data.begin() + (i * numComps);
00186   while (numComps-- > 0)
00187     *(tuple++) = static_cast<double>(*(in++));
00188 }
00189 
00190 template <typename Scalar>
00191 vtkIdType TeuchosArrayRCPDataArray<Scalar>::LookupTypedValue(ValueType value)
00192 {
00193   return this->Lookup(value, 0);
00194 }
00195 
00196 template <typename Scalar>
00197 void TeuchosArrayRCPDataArray<Scalar>::
00198 LookupTypedValue(ValueType value, vtkIdList *ids)
00199 {
00200   ids->Reset();
00201   vtkIdType index = 0;
00202   while ((index = this->Lookup(value, index)) >= 0)
00203     ids->InsertNextId(index++);
00204 }
00205 
00206 template <typename Scalar>
00207 int TeuchosArrayRCPDataArray<Scalar>::Allocate(vtkIdType sz, vtkIdType ext)
00208 {
00209   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00210   return 0;
00211 }
00212 
00213 template <typename Scalar>
00214 int TeuchosArrayRCPDataArray<Scalar>::Resize(vtkIdType numTuples)
00215 {
00216   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00217   return 0;
00218 }
00219 
00220 template <typename Scalar>
00221 void TeuchosArrayRCPDataArray<Scalar>::SetNumberOfTuples(vtkIdType number)
00222 {
00223   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00224 }
00225 
00226 template <typename Scalar>
00227 void TeuchosArrayRCPDataArray<Scalar>::
00228 SetTuple(vtkIdType i, vtkIdType j, vtkAbstractArray *source)
00229 {
00230   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00231 }
00232 
00233 template <typename Scalar>
00234 void TeuchosArrayRCPDataArray<Scalar>::SetTuple(vtkIdType i,
00235                                                 const float *source)
00236 {
00237   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00238 }
00239 
00240 template <typename Scalar>
00241 void TeuchosArrayRCPDataArray<Scalar>::SetTuple(vtkIdType i,
00242                                                 const double *source)
00243 {
00244   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00245 }
00246 
00247 template <typename Scalar>
00248 void TeuchosArrayRCPDataArray<Scalar>::
00249 InsertTuple(vtkIdType i, vtkIdType j, vtkAbstractArray *source)
00250 {
00251   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00252 }
00253 
00254 template <typename Scalar>
00255 void TeuchosArrayRCPDataArray<Scalar>::
00256 InsertTuple(vtkIdType i, const float *source)
00257 {
00258   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00259 }
00260 
00261 template <typename Scalar>
00262 void TeuchosArrayRCPDataArray<Scalar>::
00263 InsertTuple(vtkIdType i, const double *source)
00264 {
00265   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00266 }
00267 
00268 template <typename Scalar>
00269 void TeuchosArrayRCPDataArray<Scalar>::
00270 InsertTuples(vtkIdList *dstIds, vtkIdList *srcIds, vtkAbstractArray *source)
00271 {
00272   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00273 }
00274 
00275 template <typename Scalar>
00276 vtkIdType TeuchosArrayRCPDataArray<Scalar>::
00277 InsertNextTuple(vtkIdType j, vtkAbstractArray *source)
00278 {
00279   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00280   return -1;
00281 }
00282 
00283 template <typename Scalar>
00284 vtkIdType TeuchosArrayRCPDataArray<Scalar>::InsertNextTuple(const float *source)
00285 {
00286   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00287   return -1;
00288 }
00289 
00290 template <typename Scalar>
00291 vtkIdType TeuchosArrayRCPDataArray<Scalar>::
00292 InsertNextTuple(const double *source)
00293 {
00294   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00295   return -1;
00296 }
00297 
00298 template <typename Scalar>
00299 void TeuchosArrayRCPDataArray<Scalar>::DeepCopy(vtkAbstractArray *aa)
00300 {
00301   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00302 }
00303 
00304 template <typename Scalar>
00305 void TeuchosArrayRCPDataArray<Scalar>::DeepCopy(vtkDataArray *da)
00306 {
00307   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00308 }
00309 
00310 template <typename Scalar>
00311 void TeuchosArrayRCPDataArray<Scalar>::
00312 InterpolateTuple(vtkIdType i, vtkIdList *ptIndices, vtkAbstractArray *source,
00313                  double *weights)
00314 {
00315   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00316 }
00317 
00318 template <typename Scalar>
00319 void TeuchosArrayRCPDataArray<Scalar>::
00320 InterpolateTuple(vtkIdType i, vtkIdType id1, vtkAbstractArray *source1,
00321                  vtkIdType id2, vtkAbstractArray *source2, double t)
00322 {
00323   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00324 }
00325 
00326 template <typename Scalar>
00327 void TeuchosArrayRCPDataArray<Scalar>::
00328 SetVariantValue(vtkIdType idx, vtkVariant value)
00329 {
00330   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00331 }
00332 
00333 template <typename Scalar>
00334 void TeuchosArrayRCPDataArray<Scalar>::RemoveTuple(vtkIdType id)
00335 {
00336   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00337 }
00338 
00339 template <typename Scalar>
00340 void TeuchosArrayRCPDataArray<Scalar>::RemoveFirstTuple()
00341 {
00342   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00343 }
00344 
00345 template <typename Scalar>
00346 void TeuchosArrayRCPDataArray<Scalar>::RemoveLastTuple()
00347 {
00348   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00349 }
00350 
00351 template <typename Scalar>
00352 void TeuchosArrayRCPDataArray<Scalar>::
00353 SetTupleValue(vtkIdType i, const ValueType *t)
00354 {
00355   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00356 }
00357 
00358 template <typename Scalar>
00359 void TeuchosArrayRCPDataArray<Scalar>::
00360 InsertTupleValue(vtkIdType i, const ValueType *t)
00361 {
00362   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00363 }
00364 
00365 template <typename Scalar>
00366 vtkIdType TeuchosArrayRCPDataArray<Scalar>::
00367 InsertNextTupleValue(const ValueType *t)
00368 {
00369   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00370   return -1;
00371 }
00372 
00373 template <typename Scalar>
00374 void TeuchosArrayRCPDataArray<Scalar>::SetValue(vtkIdType idx, ValueType value)
00375 {
00376   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00377 }
00378 
00379 template <typename Scalar>
00380 vtkIdType TeuchosArrayRCPDataArray<Scalar>::InsertNextValue(ValueType v)
00381 {
00382   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00383   return -1;
00384 }
00385 
00386 template <typename Scalar>
00387 void TeuchosArrayRCPDataArray<Scalar>::InsertValue(vtkIdType idx, ValueType v)
00388 {
00389   TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Read-only container.");
00390 }
00391 
00392 template <typename Scalar>
00393 TeuchosArrayRCPDataArray<Scalar>::TeuchosArrayRCPDataArray()
00394   : Data(Teuchos::ArrayRCP<Scalar>()),
00395     TmpArray(NULL)
00396 {
00397 }
00398 
00399 template <typename Scalar>
00400 TeuchosArrayRCPDataArray<Scalar>::~TeuchosArrayRCPDataArray()
00401 {
00402   delete TmpArray;
00403 }
00404 
00405 template <typename Scalar>
00406 vtkIdType TeuchosArrayRCPDataArray<Scalar>::
00407 Lookup(double val, vtkIdType index)
00408 {
00409   while (index <= this->MaxId) {
00410     if (this->Data[index] == val)
00411       return index;
00412     ++index;
00413   }
00414   return -1;
00415 }
00416 
00417 } // namespace Catalyst
00418 } // namespace Albany

Generated on Wed Mar 26 2014 18:36:36 for Albany: a Trilinos-based PDE code by  doxygen 1.7.1