/[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

trunk/esys2/escript/src/Data/DataExpanded.cpp revision 126 by jgs, Fri Jul 22 03:53:08 2005 UTC trunk/escript/src/DataExpanded.cpp revision 983 by gross, Tue Feb 20 02:49:08 2007 UTC
# Line 1  Line 1 
1  // $Id$  // $Id$
2  /*  /*
3   ******************************************************************************   ************************************************************
4   *                                                                            *   *          Copyright 2006 by ACcESS MNRF                   *
5   *       COPYRIGHT  ACcESS 2004 -  All Rights Reserved                        *   *                                                          *
6   *                                                                            *   *              http://www.access.edu.au                    *
7   * This software is the property of ACcESS. No part of this code              *   *       Primary Business: Queensland, Australia            *
8   * may be copied in any form or by any means without the expressed written    *   *  Licensed under the Open Software License version 3.0    *
9   * consent of ACcESS.  Copying, use or modification of this software          *   *     http://www.opensource.org/licenses/osl-3.0.php       *
10   * by any unauthorised person is illegal unless that person has a software    *   *                                                          *
11   * license agreement with ACcESS.                                             *   ************************************************************
  *                                                                            *  
  ******************************************************************************  
12  */  */
13    
14  #include "escript/Data/DataException.h"  #include "DataExpanded.h"
15  #include "escript/Data/DataExpanded.h"  #include "DataException.h"
16  #include "escript/Data/DataConstant.h"  #include "DataConstant.h"
17  #include "escript/Data/DataTagged.h"  #include "DataTagged.h"
18  #include "escript/Data/DataArrayView.h"  #include <netcdfcpp.h>
19    
20  #include <boost/python/extract.hpp>  #include <boost/python/extract.hpp>
21    
 #include <iostream>  
   
22  using namespace std;  using namespace std;
23  using namespace boost::python;  using namespace boost::python;
24  using namespace boost;  using namespace boost;
# Line 148  DataExpanded::DataExpanded(const Functio Line 144  DataExpanded::DataExpanded(const Functio
144    : DataAbstract(what)    : DataAbstract(what)
145  {  {
146    //    //
147      // create the view of the data
148      initialise(shape,what.getNumSamples(),what.getNumDPPSample());
149      //
150    // copy the data in the correct format    // copy the data in the correct format
151    m_data.getData()=data;    m_data.getData()=data;
   //  
   // create the view of the data  
   DataArrayView tempView(m_data.getData(),shape);  
   setPointDataView(tempView);  
152  }  }
153    
154  DataExpanded::~DataExpanded()  DataExpanded::~DataExpanded()
155  {  {
156  }  }
157    
 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);  
 }  
   
158  DataAbstract*  DataAbstract*
159  DataExpanded::getSlice(const DataArrayView::RegionType& region) const  DataExpanded::getSlice(const DataArrayView::RegionType& region) const
160  {  {
# Line 219  DataExpanded::setSlice(const DataAbstrac Line 178  DataExpanded::setSlice(const DataAbstrac
178    if (getPointDataView().getRank()!=region.size()) {    if (getPointDataView().getRank()!=region.size()) {
179      throw DataException("Error - Invalid slice region.");      throw DataException("Error - Invalid slice region.");
180    }    }
181    if (tempDataExp->getPointDataView().getRank()>0 and !value->getPointDataView().checkShape(shape)) {    if (tempDataExp->getPointDataView().getRank()>0 && !value->getPointDataView().checkShape(shape)) {
182      throw DataException (value->getPointDataView().createShapeErrorMessage(      throw DataException (value->getPointDataView().createShapeErrorMessage(
183          "Error - Couldn't copy slice due to shape mismatch.",shape));          "Error - Couldn't copy slice due to shape mismatch.",shape));
184    }    }
# Line 295  string Line 254  string
254  DataExpanded::toString() const  DataExpanded::toString() const
255  {  {
256    stringstream temp;    stringstream temp;
257      FunctionSpace fs=getFunctionSpace();
258    //    //
259    // create a temporary view as the offset will be changed    // create a temporary view as the offset will be changed
260    DataArrayView tempView(getPointDataView().getData(),getPointDataView().getShape(),getPointDataView().getOffset());    DataArrayView tempView(getPointDataView().getData(),getPointDataView().getShape(),getPointDataView().getOffset());
# Line 302  DataExpanded::toString() const Line 262  DataExpanded::toString() const
262      for (int j=0;j<m_data.getNumCols();j++) {      for (int j=0;j<m_data.getNumCols();j++) {
263        tempView.setOffset(m_data.index(i,j));        tempView.setOffset(m_data.index(i,j));
264        stringstream suffix;        stringstream suffix;
265        suffix << "(" << i << "," << j << ")";        suffix << "( id: " << i << ", ref: " << fs.getReferenceIDOfSample(i) << ", pnt: " << j << ")";
266        temp << tempView.toString(suffix.str());        temp << tempView.toString(suffix.str());
267        if (!(i==(m_data.getNumRows()-1) && j==(m_data.getNumCols()-1))) {        if (!(i==(m_data.getNumRows()-1) && j==(m_data.getNumCols()-1))) {
268          temp << endl;          temp << endl;
# Line 333  DataExpanded::getLength() const Line 293  DataExpanded::getLength() const
293    return m_data.size();    return m_data.size();
294  }  }
295    
 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);  
 }  
   
296  int  int
297  DataExpanded::archiveData(ofstream& archiveFile,  DataExpanded::archiveData(ofstream& archiveFile,
298                            const DataArrayView::ValueType::size_type noValues) const                            const DataArrayView::ValueType::size_type noValues) const
# Line 414  DataExpanded::extractData(ifstream& arch Line 307  DataExpanded::extractData(ifstream& arch
307    return(m_data.extractData(archiveFile, noValues));    return(m_data.extractData(archiveFile, noValues));
308  }  }
309    
310    void
311    DataExpanded::copyToDataPoint(const int sampleNo, const int dataPointNo, const double value) {
312      //
313      // Get the number of samples and data-points per sample.
314      int numSamples = getNumSamples();
315      int numDataPointsPerSample = getNumDPPSample();
316      int dataPointRank = getPointDataView().getRank();
317      ShapeType dataPointShape = getPointDataView().getShape();
318      if (numSamples*numDataPointsPerSample > 0) {
319         //TODO: global error handling
320         if ((sampleNo >= numSamples) || (sampleNo < 0 )) {
321              throw DataException("Error - DataExpanded::copyDataPoint invalid sampleNo.");
322         }
323         if ((dataPointNo >= numDataPointsPerSample) || (dataPointNo < 0)) {
324               throw DataException("Error - DataExpanded::copyDataPoint invalid dataPointNoInSample.");
325         }
326         DataArrayView dataPointView = getDataPoint(sampleNo, dataPointNo);
327         if (dataPointRank==0) {
328             dataPointView()=value;
329         } else if (dataPointRank==1) {
330            for (int i=0; i<dataPointShape[0]; i++) {
331                dataPointView(i)=value;
332            }
333         } else if (dataPointRank==2) {
334            for (int i=0; i<dataPointShape[0]; i++) {
335               for (int j=0; j<dataPointShape[1]; j++) {
336                  dataPointView(i,j)=value;
337               }
338            }
339         } else if (dataPointRank==3) {
340            for (int i=0; i<dataPointShape[0]; i++) {
341               for (int j=0; j<dataPointShape[1]; j++) {
342                  for (int k=0; k<dataPointShape[2]; k++) {
343                     dataPointView(i,j,k)=value;
344                  }
345               }
346            }
347         } else if (dataPointRank==4) {
348             for (int i=0; i<dataPointShape[0]; i++) {
349               for (int j=0; j<dataPointShape[1]; j++) {
350                 for (int k=0; k<dataPointShape[2]; k++) {
351                   for (int l=0; l<dataPointShape[3]; l++) {
352                      dataPointView(i,j,k,l)=value;
353                   }
354                 }
355               }
356             }
357         }
358      }
359    }
360    void
361    DataExpanded::copyToDataPoint(const int sampleNo, const int dataPointNo, const boost::python::numeric::array& value) {
362      //
363      // Get the number of samples and data-points per sample.
364      int numSamples = getNumSamples();
365      int numDataPointsPerSample = getNumDPPSample();
366      int dataPointRank = getPointDataView().getRank();
367      ShapeType dataPointShape = getPointDataView().getShape();
368      //
369      // check rank:
370      if (value.getrank()!=dataPointRank)
371           throw DataException("Rank of numarray does not match Data object rank");
372      if (numSamples*numDataPointsPerSample > 0) {
373         //TODO: global error handling
374         if ((sampleNo >= numSamples) || (sampleNo < 0 )) {
375              throw DataException("Error - DataExpanded::copyDataPoint invalid sampleNo.");
376         }
377         if ((dataPointNo >= numDataPointsPerSample) || (dataPointNo < 0)) {
378               throw DataException("Error - DataExpanded::copyDataPoint invalid dataPointNoInSample.");
379         }
380         DataArrayView dataPointView = getDataPoint(sampleNo, dataPointNo);
381         if (dataPointRank==0) {
382             dataPointView()=extract<double>(value[0]);
383         } else if (dataPointRank==1) {
384            for (int i=0; i<dataPointShape[0]; i++) {
385                dataPointView(i)=extract<double>(value[i]);
386            }
387         } else if (dataPointRank==2) {
388            for (int i=0; i<dataPointShape[0]; i++) {
389               for (int j=0; j<dataPointShape[1]; j++) {
390                  dataPointView(i,j)=extract<double>(value[i][j]);
391               }
392            }
393         } else if (dataPointRank==3) {
394            for (int i=0; i<dataPointShape[0]; i++) {
395               for (int j=0; j<dataPointShape[1]; j++) {
396                  for (int k=0; k<dataPointShape[2]; k++) {
397                     dataPointView(i,j,k)=extract<double>(value[i][j][k]);
398                  }
399               }
400            }
401         } else if (dataPointRank==4) {
402             for (int i=0; i<dataPointShape[0]; i++) {
403               for (int j=0; j<dataPointShape[1]; j++) {
404                 for (int k=0; k<dataPointShape[2]; k++) {
405                   for (int l=0; l<dataPointShape[3]; l++) {
406                      dataPointView(i,j,k,l)=extract<double>(value[i][j][k][l]);
407                   }
408                 }
409               }
410             }
411         }
412      }
413    }
414  void  void
415  DataExpanded::copyAll(const boost::python::numeric::array& value) {  DataExpanded::copyAll(const boost::python::numeric::array& value) {
416    //    //
# Line 468  DataExpanded::copyAll(const boost::pytho Line 465  DataExpanded::copyAll(const boost::pytho
465      }      }
466    }    }
467  }  }
468    void
469    DataExpanded::symmetric(DataAbstract* ev)
470    {
471      int sampleNo,dataPointNo;
472      int numSamples = getNumSamples();
473      int numDataPointsPerSample = getNumDPPSample();
474      DataExpanded* temp_ev=dynamic_cast<DataExpanded*>(ev);
475      if (temp_ev==0) {
476        throw DataException("Error - DataExpanded::symmetric: casting to DataExpanded failed (propably a programming error).");
477      }
478      DataArrayView& thisView=getPointDataView();
479      DataArrayView& evView=ev->getPointDataView();
480      #pragma omp parallel for private(sampleNo,dataPointNo) schedule(static)
481      for (sampleNo = 0; sampleNo < numSamples; sampleNo++) {
482        for (dataPointNo = 0; dataPointNo < numDataPointsPerSample; dataPointNo++) {
483             DataArrayView::symmetric(thisView,getPointOffset(sampleNo,dataPointNo),
484                                        evView,ev->getPointOffset(sampleNo,dataPointNo));
485        }
486      }
487    }
488    void
489    DataExpanded::nonsymmetric(DataAbstract* ev)
490    {
491      int sampleNo,dataPointNo;
492      int numSamples = getNumSamples();
493      int numDataPointsPerSample = getNumDPPSample();
494      DataExpanded* temp_ev=dynamic_cast<DataExpanded*>(ev);
495      if (temp_ev==0) {
496        throw DataException("Error - DataExpanded::nonsymmetric: casting to DataExpanded failed (propably a programming error).");
497      }
498      DataArrayView& thisView=getPointDataView();
499      DataArrayView& evView=ev->getPointDataView();
500      #pragma omp parallel for private(sampleNo,dataPointNo) schedule(static)
501      for (sampleNo = 0; sampleNo < numSamples; sampleNo++) {
502        for (dataPointNo = 0; dataPointNo < numDataPointsPerSample; dataPointNo++) {
503             DataArrayView::nonsymmetric(thisView,getPointOffset(sampleNo,dataPointNo),
504                                        evView,ev->getPointOffset(sampleNo,dataPointNo));
505        }
506      }
507    }
508    void
509    DataExpanded::trace(DataAbstract* ev, int axis_offset)
510    {
511      int sampleNo,dataPointNo;
512      int numSamples = getNumSamples();
513      int numDataPointsPerSample = getNumDPPSample();
514      DataExpanded* temp_ev=dynamic_cast<DataExpanded*>(ev);
515      if (temp_ev==0) {
516        throw DataException("Error - DataExpanded::trace: casting to DataExpanded failed (propably a programming error).");
517      }
518      DataArrayView& thisView=getPointDataView();
519      DataArrayView& evView=ev->getPointDataView();
520      #pragma omp parallel for private(sampleNo,dataPointNo) schedule(static)
521      for (sampleNo = 0; sampleNo < numSamples; sampleNo++) {
522        for (dataPointNo = 0; dataPointNo < numDataPointsPerSample; dataPointNo++) {
523             DataArrayView::trace(thisView,getPointOffset(sampleNo,dataPointNo),
524                                        evView,ev->getPointOffset(sampleNo,dataPointNo),axis_offset);
525        }
526      }
527    }
528    
529    void
530    DataExpanded::transpose(DataAbstract* ev, int axis_offset)
531    {
532      int sampleNo,dataPointNo;
533      int numSamples = getNumSamples();
534      int numDataPointsPerSample = getNumDPPSample();
535      DataExpanded* temp_ev=dynamic_cast<DataExpanded*>(ev);
536      if (temp_ev==0) {
537        throw DataException("Error - DataExpanded::transpose: casting to DataExpanded failed (propably a programming error).");
538      }
539      DataArrayView& thisView=getPointDataView();
540      DataArrayView& evView=ev->getPointDataView();
541      #pragma omp parallel for private(sampleNo,dataPointNo) schedule(static)
542      for (sampleNo = 0; sampleNo < numSamples; sampleNo++) {
543        for (dataPointNo = 0; dataPointNo < numDataPointsPerSample; dataPointNo++) {
544             DataArrayView::transpose(thisView,getPointOffset(sampleNo,dataPointNo),
545                                        evView,ev->getPointOffset(sampleNo,dataPointNo),axis_offset);
546        }
547      }
548    }
549    
550    void
551    DataExpanded::swapaxes(DataAbstract* ev, int axis0, int axis1)
552    {
553      int sampleNo,dataPointNo;
554      int numSamples = getNumSamples();
555      int numDataPointsPerSample = getNumDPPSample();
556      DataExpanded* temp_ev=dynamic_cast<DataExpanded*>(ev);
557      if (temp_ev==0) {
558        throw DataException("Error - DataExpanded::swapaxes: casting to DataExpanded failed (propably a programming error).");
559      }
560      DataArrayView& thisView=getPointDataView();
561      DataArrayView& evView=ev->getPointDataView();
562      #pragma omp parallel for private(sampleNo,dataPointNo) schedule(static)
563      for (sampleNo = 0; sampleNo < numSamples; sampleNo++) {
564        for (dataPointNo = 0; dataPointNo < numDataPointsPerSample; dataPointNo++) {
565             DataArrayView::swapaxes(thisView,getPointOffset(sampleNo,dataPointNo),
566                                        evView,ev->getPointOffset(sampleNo,dataPointNo),axis0,axis1);
567        }
568      }
569    }
570    void
571    DataExpanded::eigenvalues(DataAbstract* ev)
572    {
573      int sampleNo,dataPointNo;
574      int numSamples = getNumSamples();
575      int numDataPointsPerSample = getNumDPPSample();
576      DataExpanded* temp_ev=dynamic_cast<DataExpanded*>(ev);
577      if (temp_ev==0) {
578        throw DataException("Error - DataExpanded::eigenvalues: casting to DataExpanded failed (propably a programming error).");
579      }
580      DataArrayView& thisView=getPointDataView();
581      DataArrayView& evView=ev->getPointDataView();
582      #pragma omp parallel for private(sampleNo,dataPointNo) schedule(static)
583      for (sampleNo = 0; sampleNo < numSamples; sampleNo++) {
584        for (dataPointNo = 0; dataPointNo < numDataPointsPerSample; dataPointNo++) {
585             DataArrayView::eigenvalues(thisView,getPointOffset(sampleNo,dataPointNo),
586                                        evView,ev->getPointOffset(sampleNo,dataPointNo));
587        }
588      }
589    }
590    void
591    DataExpanded::eigenvalues_and_eigenvectors(DataAbstract* ev,DataAbstract* V,const double tol)
592    {
593      int numSamples = getNumSamples();
594      int numDataPointsPerSample = getNumDPPSample();
595      int sampleNo,dataPointNo;
596      DataExpanded* temp_ev=dynamic_cast<DataExpanded*>(ev);
597      if (temp_ev==0) {
598        throw DataException("Error - DataExpanded::eigenvalues_and_eigenvectors: casting to DataExpanded failed (propably a programming error).");
599      }
600      DataExpanded* temp_V=dynamic_cast<DataExpanded*>(V);
601      if (temp_V==0) {
602        throw DataException("Error - DataExpanded::eigenvalues_and_eigenvectors: casting to DataExpanded failed (propably a programming error).");
603      }
604      DataArrayView& thisView=getPointDataView();
605      DataArrayView& evView=ev->getPointDataView();
606      DataArrayView& VView=V->getPointDataView();
607      #pragma omp parallel for private(sampleNo,dataPointNo) schedule(static)
608      for (sampleNo = 0; sampleNo < numSamples; sampleNo++) {
609        for (dataPointNo = 0; dataPointNo < numDataPointsPerSample; dataPointNo++) {
610             DataArrayView::eigenvalues_and_eigenvectors(thisView,getPointOffset(sampleNo,dataPointNo),
611                                                         evView,ev->getPointOffset(sampleNo,dataPointNo),
612                                                         VView,V->getPointOffset(sampleNo,dataPointNo),
613                                                         tol);
614        }
615      }
616    }
617    
618    void
619    DataExpanded::dump(const std::string fileName) const
620    {
621       #ifdef PASO_MPI
622       throw DataException("Error - DataExpanded:: dump is not implemented for MPI yet.");
623       #endif
624       const int ldims=2+DataArrayView::maxRank;
625       const NcDim* ncdims[ldims];
626       NcVar *var, *ids;
627       int rank = getPointDataView().getRank();
628       int type=  getFunctionSpace().getTypeCode();
629       int ndims =0;
630       long dims[ldims];
631       DataArrayView::ShapeType shape = getPointDataView().getShape();
632    
633       // netCDF error handler
634       NcError err(NcError::verbose_nonfatal);
635       // Create the file.
636       NcFile dataFile(fileName.c_str(), NcFile::Replace);
637       // check if writing was successful
638       if (!dataFile.is_valid())
639            throw DataException("Error - DataExpanded:: opening of netCDF file for output failed.");
640       if (!dataFile.add_att("type","expanded") )
641            throw DataException("Error - DataExpanded:: appending data type to netCDF file failed.");
642       if (!dataFile.add_att("rank",rank) )
643            throw DataException("Error - DataExpanded:: appending rank attribute to netCDF file failed.");
644       if (!dataFile.add_att("function_space_type",type))
645            throw DataException("Error - DataExpanded:: appending function space attribute to netCDF file failed.");
646       ndims=rank+2;
647       if ( rank >0 ) {
648           dims[0]=shape[0];
649           if (! (ncdims[0] = dataFile.add_dim("d0",shape[0])) )
650                throw DataException("Error - DataExpanded:: appending ncdimsion 0 to netCDF file failed.");
651       }
652       if ( rank >1 ) {
653           dims[1]=shape[1];
654           if (! (ncdims[1] = dataFile.add_dim("d1",shape[1])) )
655                throw DataException("Error - DataExpanded:: appending ncdimsion 1 to netCDF file failed.");
656       }
657       if ( rank >2 ) {
658           dims[2]=shape[2];
659           if (! (ncdims[2] = dataFile.add_dim("d2", shape[2])) )
660                throw DataException("Error - DataExpanded:: appending ncdimsion 2 to netCDF file failed.");
661       }
662       if ( rank >3 ) {
663           dims[3]=shape[3];
664           if (! (ncdims[3] = dataFile.add_dim("d3", shape[3])) )
665                throw DataException("Error - DataExpanded:: appending ncdimsion 3 to netCDF file failed.");
666       }
667       dims[rank]=getFunctionSpace().getNumDataPointsPerSample();
668       if (! (ncdims[rank] = dataFile.add_dim("num_data_points_per_sample", dims[rank])) )
669                throw DataException("Error - DataExpanded:: appending num_data_points_per_sample to netCDF file failed.");
670       dims[rank+1]=getFunctionSpace().getNumSamples();
671       if (! (ncdims[rank+1] = dataFile.add_dim("num_samples", dims[rank+1])) )
672                throw DataException("Error - DataExpanded:: appending num_sample to netCDF file failed.");
673    
674       if (! ( ids = dataFile.add_var("id", ncInt, ncdims[rank+1])) )
675            throw DataException("Error - DataExpanded:: appending reference id to netCDF file failed.");
676       const int* ids_p=getFunctionSpace().borrowSampleReferenceIDs();
677       if (! (ids->put(ids_p,dims[rank+1])) )
678            throw DataException("Error - DataExpanded:: copy reference id  to netCDF buffer failed.");
679    
680       if (! ( var = dataFile.add_var("data", ncDouble, ndims, ncdims)) )
681            throw DataException("Error - DataExpanded:: appending variable to netCDF file failed.");
682       if (! (var->put(&m_data[0],dims)) )
683            throw DataException("Error - DataExpanded:: copy data to netCDF buffer failed.");
684    }
685  }  // end of namespace  }  // end of namespace

Legend:
Removed from v.126  
changed lines
  Added in v.983

  ViewVC Help
Powered by ViewVC 1.1.26