/[escript]/trunk/escript/src/DataExpanded.cpp
ViewVC logotype

Diff of /trunk/escript/src/DataExpanded.cpp

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

revision 615 by elspeth, Wed Mar 22 02:12:00 2006 UTC revision 1118 by gross, Tue Apr 24 08:55:04 2007 UTC
# Line 15  Line 15 
15  #include "DataException.h"  #include "DataException.h"
16  #include "DataConstant.h"  #include "DataConstant.h"
17  #include "DataTagged.h"  #include "DataTagged.h"
18    #ifdef USE_NETCDF
19    #include <netcdfcpp.h>
20    #endif
21    
22  #include <boost/python/extract.hpp>  #include <boost/python/extract.hpp>
23    
# Line 154  DataExpanded::~DataExpanded() Line 157  DataExpanded::~DataExpanded()
157  {  {
158  }  }
159    
 void  
 DataExpanded::reshapeDataPoint(const DataArrayView::ShapeType& shape)  
 {  
   if (getPointDataView().getRank()!=0) {  
     stringstream temp;  
     temp << "Error - Can only reshape Data with data points of rank 0. "  
          << "This Data has data points with rank: "  
          << getPointDataView().getRank();  
     throw DataException(temp.str());  
   }  
   //  
   // create the new DataBlocks2D data array, and a corresponding DataArrayView  
   DataBlocks2D newData(getNumSamples(),getNumDPPSample(),DataArrayView::noValues(shape));  
   DataArrayView newView(newData.getData(),shape);  
   //  
   // Copy the original data to every value for the new shape  
   int i,j;  
   int nRows=m_data.getNumRows();  
   int nCols=m_data.getNumCols();  
   #pragma omp parallel for private(i,j) schedule(static)  
   for (i=0;i<nRows;i++) {  
     for (j=0;j<nCols;j++) {  
       // NOTE: An exception may be thown from this call if  
       // DOASSERT is on which of course will play  
       // havoc with the omp threads. Run single threaded  
       // if using DOASSERT.  
       newView.copy(newData.index(i,j),m_data(i,j));  
     }  
   }  
   // swap the new data array for the original  
   m_data.Swap(newData);  
   // set the corresponding DataArrayView  
   DataArrayView temp(m_data.getData(),shape);  
   setPointDataView(temp);  
 }  
   
160  DataAbstract*  DataAbstract*
161  DataExpanded::getSlice(const DataArrayView::RegionType& region) const  DataExpanded::getSlice(const DataArrayView::RegionType& region) const
162  {  {
# Line 213  DataExpanded::setSlice(const DataAbstrac Line 180  DataExpanded::setSlice(const DataAbstrac
180    if (getPointDataView().getRank()!=region.size()) {    if (getPointDataView().getRank()!=region.size()) {
181      throw DataException("Error - Invalid slice region.");      throw DataException("Error - Invalid slice region.");
182    }    }
183    if (tempDataExp->getPointDataView().getRank()>0 and !value->getPointDataView().checkShape(shape)) {    if (tempDataExp->getPointDataView().getRank()>0 && !value->getPointDataView().checkShape(shape)) {
184      throw DataException (value->getPointDataView().createShapeErrorMessage(      throw DataException (value->getPointDataView().createShapeErrorMessage(
185          "Error - Couldn't copy slice due to shape mismatch.",shape));          "Error - Couldn't copy slice due to shape mismatch.",shape));
186    }    }
# Line 289  string Line 256  string
256  DataExpanded::toString() const  DataExpanded::toString() const
257  {  {
258    stringstream temp;    stringstream temp;
259      FunctionSpace fs=getFunctionSpace();
260    //    //
261    // create a temporary view as the offset will be changed    // create a temporary view as the offset will be changed
262    DataArrayView tempView(getPointDataView().getData(),getPointDataView().getShape(),getPointDataView().getOffset());    DataArrayView tempView(getPointDataView().getData(),getPointDataView().getShape(),getPointDataView().getOffset());
# Line 296  DataExpanded::toString() const Line 264  DataExpanded::toString() const
264      for (int j=0;j<m_data.getNumCols();j++) {      for (int j=0;j<m_data.getNumCols();j++) {
265        tempView.setOffset(m_data.index(i,j));        tempView.setOffset(m_data.index(i,j));
266        stringstream suffix;        stringstream suffix;
267        suffix << "(" << i << "," << j << ")";        suffix << "( id: " << i << ", ref: " << fs.getReferenceIDOfSample(i) << ", pnt: " << j << ")";
268        temp << tempView.toString(suffix.str());        temp << tempView.toString(suffix.str());
269        if (!(i==(m_data.getNumRows()-1) && j==(m_data.getNumCols()-1))) {        if (!(i==(m_data.getNumRows()-1) && j==(m_data.getNumCols()-1))) {
270          temp << endl;          temp << endl;
# Line 327  DataExpanded::getLength() const Line 295  DataExpanded::getLength() const
295    return m_data.size();    return m_data.size();
296  }  }
297    
 void  
 DataExpanded::setRefValue(int ref,  
                           const DataArray& value)  
 {  
   //  
   // Get the number of samples and data-points per sample.  
   int numSamples = getNumSamples();  
   int numDPPSample = getNumDPPSample();  
   
   //  
   // Determine the sample number which corresponds to this reference number.  
   int sampleNo = -1;  
   int tempRef = -1;  
   for (int n=0; n<numSamples; n++) {  
     tempRef = getFunctionSpace().getReferenceNoFromSampleNo(n);  
     if (tempRef == ref) {  
       sampleNo = n;  
       break;  
     }  
   }  
   if (sampleNo == -1) {  
     throw DataException("DataExpanded::setRefValue error: invalid ref number supplied.");  
   }  
   
   for (int n=0; n<numDPPSample; n++) {  
     //  
     // Get *each* data-point in the sample in turn.  
     DataArrayView pointView = getDataPoint(sampleNo, n);  
     //  
     // Assign the values in the DataArray to this data-point.  
     pointView.copy(value.getView());  
   }  
 }  
   
 void  
 DataExpanded::getRefValue(int ref,  
                           DataArray& value)  
 {  
   //  
   // Get the number of samples and data-points per sample.  
   int numSamples = getNumSamples();  
   int numDPPSample = getNumDPPSample();  
   
   //  
   // Determine the sample number which corresponds to this reference number  
   int sampleNo = -1;  
   int tempRef = -1;  
   for (int n=0; n<numSamples; n++) {  
     tempRef = getFunctionSpace().getReferenceNoFromSampleNo(n);  
     if (tempRef == ref) {  
       sampleNo = n;  
       break;  
     }  
   }  
   if (sampleNo == -1) {  
     throw DataException("DataExpanded::getRefValue error: invalid ref number supplied.");  
   }  
   
   //  
   // Get the *first* data-point associated with this sample number.  
   DataArrayView pointView = getDataPoint(sampleNo, 0);  
   
   //  
   // Load the values from this data-point into the DataArray  
   value.getView().copy(pointView);  
 }  
   
298  int  int
299  DataExpanded::archiveData(ofstream& archiveFile,  DataExpanded::archiveData(ofstream& archiveFile,
300                            const DataArrayView::ValueType::size_type noValues) const                            const DataArrayView::ValueType::size_type noValues) const
# Line 408  DataExpanded::extractData(ifstream& arch Line 309  DataExpanded::extractData(ifstream& arch
309    return(m_data.extractData(archiveFile, noValues));    return(m_data.extractData(archiveFile, noValues));
310  }  }
311    
312    void
313    DataExpanded::copyToDataPoint(const int sampleNo, const int dataPointNo, const double value) {
314      //
315      // Get the number of samples and data-points per sample.
316      int numSamples = getNumSamples();
317      int numDataPointsPerSample = getNumDPPSample();
318      int dataPointRank = getPointDataView().getRank();
319      ShapeType dataPointShape = getPointDataView().getShape();
320      if (numSamples*numDataPointsPerSample > 0) {
321         //TODO: global error handling
322         if ((sampleNo >= numSamples) || (sampleNo < 0 )) {
323              throw DataException("Error - DataExpanded::copyDataPoint invalid sampleNo.");
324         }
325         if ((dataPointNo >= numDataPointsPerSample) || (dataPointNo < 0)) {
326               throw DataException("Error - DataExpanded::copyDataPoint invalid dataPointNoInSample.");
327         }
328         DataArrayView dataPointView = getDataPoint(sampleNo, dataPointNo);
329         if (dataPointRank==0) {
330             dataPointView()=value;
331         } else if (dataPointRank==1) {
332            for (int i=0; i<dataPointShape[0]; i++) {
333                dataPointView(i)=value;
334            }
335         } else if (dataPointRank==2) {
336            for (int i=0; i<dataPointShape[0]; i++) {
337               for (int j=0; j<dataPointShape[1]; j++) {
338                  dataPointView(i,j)=value;
339               }
340            }
341         } else if (dataPointRank==3) {
342            for (int i=0; i<dataPointShape[0]; i++) {
343               for (int j=0; j<dataPointShape[1]; j++) {
344                  for (int k=0; k<dataPointShape[2]; k++) {
345                     dataPointView(i,j,k)=value;
346                  }
347               }
348            }
349         } else if (dataPointRank==4) {
350             for (int i=0; i<dataPointShape[0]; i++) {
351               for (int j=0; j<dataPointShape[1]; j++) {
352                 for (int k=0; k<dataPointShape[2]; k++) {
353                   for (int l=0; l<dataPointShape[3]; l++) {
354                      dataPointView(i,j,k,l)=value;
355                   }
356                 }
357               }
358             }
359         }
360      }
361    }
362    void
363    DataExpanded::copyToDataPoint(const int sampleNo, const int dataPointNo, const boost::python::numeric::array& value) {
364      //
365      // Get the number of samples and data-points per sample.
366      int numSamples = getNumSamples();
367      int numDataPointsPerSample = getNumDPPSample();
368      int dataPointRank = getPointDataView().getRank();
369      ShapeType dataPointShape = getPointDataView().getShape();
370      //
371      // check rank:
372      if (value.getrank()!=dataPointRank)
373           throw DataException("Rank of numarray does not match Data object rank");
374      if (numSamples*numDataPointsPerSample > 0) {
375         //TODO: global error handling
376         if ((sampleNo >= numSamples) || (sampleNo < 0 )) {
377              throw DataException("Error - DataExpanded::copyDataPoint invalid sampleNo.");
378         }
379         if ((dataPointNo >= numDataPointsPerSample) || (dataPointNo < 0)) {
380               throw DataException("Error - DataExpanded::copyDataPoint invalid dataPointNoInSample.");
381         }
382         DataArrayView dataPointView = getDataPoint(sampleNo, dataPointNo);
383         if (dataPointRank==0) {
384             dataPointView()=extract<double>(value[0]);
385         } else if (dataPointRank==1) {
386            for (int i=0; i<dataPointShape[0]; i++) {
387                dataPointView(i)=extract<double>(value[i]);
388            }
389         } else if (dataPointRank==2) {
390            for (int i=0; i<dataPointShape[0]; i++) {
391               for (int j=0; j<dataPointShape[1]; j++) {
392                  dataPointView(i,j)=extract<double>(value[i][j]);
393               }
394            }
395         } else if (dataPointRank==3) {
396            for (int i=0; i<dataPointShape[0]; i++) {
397               for (int j=0; j<dataPointShape[1]; j++) {
398                  for (int k=0; k<dataPointShape[2]; k++) {
399                     dataPointView(i,j,k)=extract<double>(value[i][j][k]);
400                  }
401               }
402            }
403         } else if (dataPointRank==4) {
404             for (int i=0; i<dataPointShape[0]; i++) {
405               for (int j=0; j<dataPointShape[1]; j++) {
406                 for (int k=0; k<dataPointShape[2]; k++) {
407                   for (int l=0; l<dataPointShape[3]; l++) {
408                      dataPointView(i,j,k,l)=extract<double>(value[i][j][k][l]);
409                   }
410                 }
411               }
412             }
413         }
414      }
415    }
416  void  void
417  DataExpanded::copyAll(const boost::python::numeric::array& value) {  DataExpanded::copyAll(const boost::python::numeric::array& value) {
418    //    //
# Line 463  DataExpanded::copyAll(const boost::pytho Line 468  DataExpanded::copyAll(const boost::pytho
468    }    }
469  }  }
470  void  void
471    DataExpanded::symmetric(DataAbstract* ev)
472    {
473      int sampleNo,dataPointNo;
474      int numSamples = getNumSamples();
475      int numDataPointsPerSample = getNumDPPSample();
476      DataExpanded* temp_ev=dynamic_cast<DataExpanded*>(ev);
477      if (temp_ev==0) {
478        throw DataException("Error - DataExpanded::symmetric: casting to DataExpanded failed (propably a programming error).");
479      }
480      DataArrayView& thisView=getPointDataView();
481      DataArrayView& evView=ev->getPointDataView();
482      #pragma omp parallel for private(sampleNo,dataPointNo) schedule(static)
483      for (sampleNo = 0; sampleNo < numSamples; sampleNo++) {
484        for (dataPointNo = 0; dataPointNo < numDataPointsPerSample; dataPointNo++) {
485             DataArrayView::symmetric(thisView,getPointOffset(sampleNo,dataPointNo),
486                                        evView,ev->getPointOffset(sampleNo,dataPointNo));
487        }
488      }
489    }
490    void
491    DataExpanded::nonsymmetric(DataAbstract* ev)
492    {
493      int sampleNo,dataPointNo;
494      int numSamples = getNumSamples();
495      int numDataPointsPerSample = getNumDPPSample();
496      DataExpanded* temp_ev=dynamic_cast<DataExpanded*>(ev);
497      if (temp_ev==0) {
498        throw DataException("Error - DataExpanded::nonsymmetric: casting to DataExpanded failed (propably a programming error).");
499      }
500      DataArrayView& thisView=getPointDataView();
501      DataArrayView& evView=ev->getPointDataView();
502      #pragma omp parallel for private(sampleNo,dataPointNo) schedule(static)
503      for (sampleNo = 0; sampleNo < numSamples; sampleNo++) {
504        for (dataPointNo = 0; dataPointNo < numDataPointsPerSample; dataPointNo++) {
505             DataArrayView::nonsymmetric(thisView,getPointOffset(sampleNo,dataPointNo),
506                                        evView,ev->getPointOffset(sampleNo,dataPointNo));
507        }
508      }
509    }
510    void
511    DataExpanded::trace(DataAbstract* ev, int axis_offset)
512    {
513      int sampleNo,dataPointNo;
514      int numSamples = getNumSamples();
515      int numDataPointsPerSample = getNumDPPSample();
516      DataExpanded* temp_ev=dynamic_cast<DataExpanded*>(ev);
517      if (temp_ev==0) {
518        throw DataException("Error - DataExpanded::trace: casting to DataExpanded failed (propably a programming error).");
519      }
520      DataArrayView& thisView=getPointDataView();
521      DataArrayView& evView=ev->getPointDataView();
522      #pragma omp parallel for private(sampleNo,dataPointNo) schedule(static)
523      for (sampleNo = 0; sampleNo < numSamples; sampleNo++) {
524        for (dataPointNo = 0; dataPointNo < numDataPointsPerSample; dataPointNo++) {
525             DataArrayView::trace(thisView,getPointOffset(sampleNo,dataPointNo),
526                                        evView,ev->getPointOffset(sampleNo,dataPointNo),axis_offset);
527        }
528      }
529    }
530    
531    void
532    DataExpanded::transpose(DataAbstract* ev, int axis_offset)
533    {
534      int sampleNo,dataPointNo;
535      int numSamples = getNumSamples();
536      int numDataPointsPerSample = getNumDPPSample();
537      DataExpanded* temp_ev=dynamic_cast<DataExpanded*>(ev);
538      if (temp_ev==0) {
539        throw DataException("Error - DataExpanded::transpose: casting to DataExpanded failed (propably a programming error).");
540      }
541      DataArrayView& thisView=getPointDataView();
542      DataArrayView& evView=ev->getPointDataView();
543      #pragma omp parallel for private(sampleNo,dataPointNo) schedule(static)
544      for (sampleNo = 0; sampleNo < numSamples; sampleNo++) {
545        for (dataPointNo = 0; dataPointNo < numDataPointsPerSample; dataPointNo++) {
546             DataArrayView::transpose(thisView,getPointOffset(sampleNo,dataPointNo),
547                                        evView,ev->getPointOffset(sampleNo,dataPointNo),axis_offset);
548        }
549      }
550    }
551    
552    void
553    DataExpanded::swapaxes(DataAbstract* ev, int axis0, int axis1)
554    {
555      int sampleNo,dataPointNo;
556      int numSamples = getNumSamples();
557      int numDataPointsPerSample = getNumDPPSample();
558      DataExpanded* temp_ev=dynamic_cast<DataExpanded*>(ev);
559      if (temp_ev==0) {
560        throw DataException("Error - DataExpanded::swapaxes: casting to DataExpanded failed (propably a programming error).");
561      }
562      DataArrayView& thisView=getPointDataView();
563      DataArrayView& evView=ev->getPointDataView();
564      #pragma omp parallel for private(sampleNo,dataPointNo) schedule(static)
565      for (sampleNo = 0; sampleNo < numSamples; sampleNo++) {
566        for (dataPointNo = 0; dataPointNo < numDataPointsPerSample; dataPointNo++) {
567             DataArrayView::swapaxes(thisView,getPointOffset(sampleNo,dataPointNo),
568                                        evView,ev->getPointOffset(sampleNo,dataPointNo),axis0,axis1);
569        }
570      }
571    }
572    void
573  DataExpanded::eigenvalues(DataAbstract* ev)  DataExpanded::eigenvalues(DataAbstract* ev)
574  {  {
575    int sampleNo,dataPointNo;    int sampleNo,dataPointNo;
# Line 510  DataExpanded::eigenvalues_and_eigenvecto Line 617  DataExpanded::eigenvalues_and_eigenvecto
617    }    }
618  }  }
619    
620    void
621    DataExpanded::setToZero(){
622      int numSamples = getNumSamples();
623      int numDataPointsPerSample = getNumDPPSample();
624      DataArrayView& thisView=getPointDataView();
625      DataArrayView::ValueType::size_type n = thisView.noValues();
626      double* p;
627      int  sampleNo,dataPointNo, i;
628      #pragma omp parallel for private(sampleNo,dataPointNo,p,i) schedule(static)
629      for (sampleNo = 0; sampleNo < numSamples; sampleNo++) {
630        for (dataPointNo = 0; dataPointNo < numDataPointsPerSample; dataPointNo++) {
631            p=&(m_data[getPointOffset(sampleNo,dataPointNo)]);
632            for (int i=0; i<n ;++i) p[i]=0.;
633        }
634      }
635    }
636    
637    
638    void
639    DataExpanded::dump(const std::string fileName) const
640    {
641       #ifdef PASO_MPI
642       throw DataException("Error - DataExpanded:: dump is not implemented for MPI yet.");
643       #endif
644       #ifdef USE_NETCDF
645       const int ldims=2+DataArrayView::maxRank;
646       const NcDim* ncdims[ldims];
647       NcVar *var, *ids;
648       int rank = getPointDataView().getRank();
649       int type=  getFunctionSpace().getTypeCode();
650       int ndims =0;
651       long dims[ldims];
652       DataArrayView::ShapeType shape = getPointDataView().getShape();
653    
654       // netCDF error handler
655       NcError err(NcError::verbose_nonfatal);
656       // Create the file.
657       NcFile dataFile(fileName.c_str(), NcFile::Replace);
658       // check if writing was successful
659       if (!dataFile.is_valid())
660            throw DataException("Error - DataExpanded:: opening of netCDF file for output failed.");
661       if (!dataFile.add_att("type","expanded") )
662            throw DataException("Error - DataExpanded:: appending data type to netCDF file failed.");
663       if (!dataFile.add_att("rank",rank) )
664            throw DataException("Error - DataExpanded:: appending rank attribute to netCDF file failed.");
665       if (!dataFile.add_att("function_space_type",type))
666            throw DataException("Error - DataExpanded:: appending function space attribute to netCDF file failed.");
667       ndims=rank+2;
668       if ( rank >0 ) {
669           dims[0]=shape[0];
670           if (! (ncdims[0] = dataFile.add_dim("d0",shape[0])) )
671                throw DataException("Error - DataExpanded:: appending ncdimsion 0 to netCDF file failed.");
672       }
673       if ( rank >1 ) {
674           dims[1]=shape[1];
675           if (! (ncdims[1] = dataFile.add_dim("d1",shape[1])) )
676                throw DataException("Error - DataExpanded:: appending ncdimsion 1 to netCDF file failed.");
677       }
678       if ( rank >2 ) {
679           dims[2]=shape[2];
680           if (! (ncdims[2] = dataFile.add_dim("d2", shape[2])) )
681                throw DataException("Error - DataExpanded:: appending ncdimsion 2 to netCDF file failed.");
682       }
683       if ( rank >3 ) {
684           dims[3]=shape[3];
685           if (! (ncdims[3] = dataFile.add_dim("d3", shape[3])) )
686                throw DataException("Error - DataExpanded:: appending ncdimsion 3 to netCDF file failed.");
687       }
688       dims[rank]=getFunctionSpace().getNumDataPointsPerSample();
689       if (! (ncdims[rank] = dataFile.add_dim("num_data_points_per_sample", dims[rank])) )
690                throw DataException("Error - DataExpanded:: appending num_data_points_per_sample to netCDF file failed.");
691       dims[rank+1]=getFunctionSpace().getNumSamples();
692       if (! (ncdims[rank+1] = dataFile.add_dim("num_samples", dims[rank+1])) )
693                throw DataException("Error - DataExpanded:: appending num_sample to netCDF file failed.");
694    
695       if (! ( ids = dataFile.add_var("id", ncInt, ncdims[rank+1])) )
696            throw DataException("Error - DataExpanded:: appending reference id to netCDF file failed.");
697       const int* ids_p=getFunctionSpace().borrowSampleReferenceIDs();
698       if (! (ids->put(ids_p,dims[rank+1])) )
699            throw DataException("Error - DataExpanded:: copy reference id  to netCDF buffer failed.");
700    
701       if (! ( var = dataFile.add_var("data", ncDouble, ndims, ncdims)) )
702            throw DataException("Error - DataExpanded:: appending variable to netCDF file failed.");
703       if (! (var->put(&m_data[0],dims)) )
704            throw DataException("Error - DataExpanded:: copy data to netCDF buffer failed.");
705       #else
706       throw DataException("Error - DataExpanded:: dump is not configured with netCDF. Please contact your installation manager.");
707       #endif
708    }
709  }  // end of namespace  }  // end of namespace

Legend:
Removed from v.615  
changed lines
  Added in v.1118

  ViewVC Help
Powered by ViewVC 1.1.26