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

Annotation of /trunk/escript/src/Data.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2246 - (hide annotations)
Thu Feb 5 06:04:55 2009 UTC (11 years, 2 months ago) by jfenwick
File size: 85997 byte(s)
Made some changes to the way copyWithMask works.

1 jgs 480
2 ksteube 1312 /*******************************************************
3 ksteube 1811 *
4     * Copyright (c) 2003-2008 by University of Queensland
5     * Earth Systems Science Computational Center (ESSCC)
6     * http://www.uq.edu.au/esscc
7     *
8     * Primary Business: Queensland, Australia
9     * Licensed under the Open Software License version 3.0
10     * http://www.opensource.org/licenses/osl-3.0.php
11     *
12     *******************************************************/
13 ksteube 1312
14 ksteube 1811
15 jgs 474 #include "Data.h"
16 jgs 94
17 jgs 480 #include "DataExpanded.h"
18     #include "DataConstant.h"
19     #include "DataTagged.h"
20     #include "DataEmpty.h"
21 jfenwick 2005 #include "DataLazy.h"
22 jgs 480 #include "FunctionSpaceFactory.h"
23     #include "AbstractContinuousDomain.h"
24     #include "UnaryFuncs.h"
25 jfenwick 1796 #include "FunctionSpaceException.h"
26 jfenwick 1897 #include "EscriptParams.h"
27 jfenwick 1796
28 ksteube 1312 extern "C" {
29 phornby 2078 #include "esysUtils/blocktimer.h"
30 ksteube 1312 }
31 jgs 480
32 jgs 119 #include <fstream>
33 jgs 94 #include <algorithm>
34     #include <vector>
35     #include <functional>
36 jfenwick 2086 #include <sstream> // so we can throw messages about ranks
37 jgs 94
38 jgs 153 #include <boost/python/dict.hpp>
39 jgs 94 #include <boost/python/extract.hpp>
40     #include <boost/python/long.hpp>
41    
42     using namespace std;
43     using namespace boost::python;
44     using namespace boost;
45     using namespace escript;
46    
47 jfenwick 2005 // ensure the current object is not a DataLazy
48     // The idea was that we could add an optional warning whenever a resolve is forced
49     #define FORCERESOLVE if (isLazy()) {resolve();}
50 jfenwick 2146 #define AUTOLAZYON escriptParams.getAUTOLAZY()
51     #define MAKELAZYOP(X) if (isLazy() || (AUTOLAZYON && m_data->isExpanded())) \
52     {\
53     DataLazy* c=new DataLazy(borrowDataPtr(),X);\
54     return Data(c);\
55     }
56     #define MAKELAZYOPOFF(X,Y) if (isLazy() || (AUTOLAZYON && m_data->isExpanded())) \
57     {\
58     DataLazy* c=new DataLazy(borrowDataPtr(),X,Y);\
59     return Data(c);\
60     }
61 jfenwick 2005
62 jfenwick 2146 #define MAKELAZYBINSELF(R,X) if (isLazy() || R.isLazy() || (AUTOLAZYON && (isExpanded() || R.isExpanded()))) \
63     {\
64     DataLazy* c=new DataLazy(m_data,R.borrowDataPtr(),X);\
65     m_data=c->getPtr();\
66     return (*this);\
67     }
68    
69     // like the above but returns a new data rather than *this
70     #define MAKELAZYBIN(R,X) if (isLazy() || R.isLazy() || (AUTOLAZYON && (isExpanded() || R.isExpanded()))) \
71     {\
72     DataLazy* c=new DataLazy(m_data,R.borrowDataPtr(),X);\
73     return Data(c);\
74     }
75    
76     #define MAKELAZYBIN2(L,R,X) if (L.isLazy() || R.isLazy() || (AUTOLAZYON && (L.isExpanded() || R.isExpanded()))) \
77     {\
78     DataLazy* c=new DataLazy(L.borrowDataPtr(),R.borrowDataPtr(),X);\
79     return Data(c);\
80     }
81    
82 jgs 94 Data::Data()
83     {
84     //
85     // Default data is type DataEmpty
86     DataAbstract* temp=new DataEmpty();
87 jfenwick 1872 m_data=temp->getPtr();
88 gross 783 m_protected=false;
89 jgs 94 }
90    
91     Data::Data(double value,
92     const tuple& shape,
93     const FunctionSpace& what,
94     bool expanded)
95     {
96 jfenwick 1796 DataTypes::ShapeType dataPointShape;
97 jgs 94 for (int i = 0; i < shape.attr("__len__")(); ++i) {
98     dataPointShape.push_back(extract<const int>(shape[i]));
99     }
100 matt 1319
101 jfenwick 1796 int len = DataTypes::noValues(dataPointShape);
102 matt 1319 DataVector temp_data(len,value,len);
103 jfenwick 1796 initialise(temp_data, dataPointShape, what, expanded);
104 gross 783 m_protected=false;
105 jgs 94 }
106    
107     Data::Data(double value,
108 jfenwick 1796 const DataTypes::ShapeType& dataPointShape,
109 jgs 94 const FunctionSpace& what,
110     bool expanded)
111     {
112 jfenwick 1796 int len = DataTypes::noValues(dataPointShape);
113 matt 1319
114     DataVector temp_data(len,value,len);
115 jfenwick 1796 // DataArrayView temp_dataView(temp_data, dataPointShape);
116 matt 1319
117 jfenwick 1796 // initialise(temp_dataView, what, expanded);
118     initialise(temp_data, dataPointShape, what, expanded);
119 matt 1319
120 gross 783 m_protected=false;
121 jgs 94 }
122    
123 jgs 102 Data::Data(const Data& inData)
124 jgs 94 {
125 jgs 102 m_data=inData.m_data;
126 gross 783 m_protected=inData.isProtected();
127 jgs 94 }
128    
129 jfenwick 1796
130 jgs 94 Data::Data(const Data& inData,
131 jfenwick 1796 const DataTypes::RegionType& region)
132 jgs 94 {
133 jfenwick 2005 DataAbstract_ptr dat=inData.m_data;
134     if (inData.isLazy())
135     {
136     dat=inData.m_data->resolve();
137     }
138     else
139     {
140     dat=inData.m_data;
141     }
142 jgs 94 //
143 jgs 102 // Create Data which is a slice of another Data
144 jfenwick 2005 DataAbstract* tmp = dat->getSlice(region);
145 jfenwick 1872 m_data=DataAbstract_ptr(tmp);
146 gross 783 m_protected=false;
147 jgs 94 }
148    
149     Data::Data(const Data& inData,
150     const FunctionSpace& functionspace)
151     {
152 jfenwick 1803 if (inData.isEmpty())
153     {
154     throw DataException("Error - will not interpolate for instances of DataEmpty.");
155     }
156 jgs 94 if (inData.getFunctionSpace()==functionspace) {
157     m_data=inData.m_data;
158 jfenwick 2005 }
159     else
160     {
161    
162     if (inData.isConstant()) { // for a constant function, we just need to use the new function space
163     if (!inData.probeInterpolation(functionspace))
164     { // Even though this is constant, we still need to check whether interpolation is allowed
165 jfenwick 2087 throw FunctionSpaceException("Cannot interpolate across to the domain of the specified FunctionSpace. (DataConstant)");
166 jfenwick 2005 }
167     // if the data is not lazy, this will just be a cast to DataReady
168     DataReady_ptr dr=inData.m_data->resolve();
169     DataConstant* dc=new DataConstant(functionspace,inData.m_data->getShape(),dr->getVector());
170     m_data=DataAbstract_ptr(dc);
171 jgs 94 } else {
172 jfenwick 2005 Data tmp(0,inData.getDataPointShape(),functionspace,true);
173     // Note: Must use a reference or pointer to a derived object
174     // in order to get polymorphic behaviour. Shouldn't really
175     // be able to create an instance of AbstractDomain but that was done
176     // as a boost:python work around which may no longer be required.
177     /*const AbstractDomain& inDataDomain=inData.getDomain();*/
178     const_Domain_ptr inDataDomain=inData.getDomain();
179     if (inDataDomain==functionspace.getDomain()) {
180     inDataDomain->interpolateOnDomain(tmp,inData);
181     } else {
182     inDataDomain->interpolateACross(tmp,inData);
183     }
184     m_data=tmp.m_data;
185 jgs 94 }
186     }
187 gross 783 m_protected=false;
188 jgs 94 }
189    
190 jfenwick 1796 Data::Data(DataAbstract* underlyingdata)
191 jgs 94 {
192 jfenwick 1872 // m_data=shared_ptr<DataAbstract>(underlyingdata);
193     m_data=underlyingdata->getPtr();
194 jfenwick 1796 m_protected=false;
195 jgs 94 }
196    
197 jfenwick 2005 Data::Data(DataAbstract_ptr underlyingdata)
198     {
199     m_data=underlyingdata;
200     m_protected=false;
201     }
202    
203    
204 jgs 94 Data::Data(const numeric::array& value,
205     const FunctionSpace& what,
206     bool expanded)
207     {
208     initialise(value,what,expanded);
209 gross 783 m_protected=false;
210 jgs 94 }
211 jfenwick 1796
212     Data::Data(const DataTypes::ValueType& value,
213     const DataTypes::ShapeType& shape,
214     const FunctionSpace& what,
215     bool expanded)
216     {
217     initialise(value,shape,what,expanded);
218     m_protected=false;
219 jgs 94 }
220    
221 jfenwick 1796
222 jgs 94 Data::Data(const object& value,
223     const FunctionSpace& what,
224     bool expanded)
225     {
226     numeric::array asNumArray(value);
227     initialise(asNumArray,what,expanded);
228 gross 783 m_protected=false;
229 jgs 94 }
230    
231 matt 1319
232 jgs 94 Data::Data(const object& value,
233     const Data& other)
234     {
235 matt 1319 numeric::array asNumArray(value);
236    
237     // extract the shape of the numarray
238 jfenwick 1796 DataTypes::ShapeType tempShape=DataTypes::shapeFromNumArray(asNumArray);
239 jfenwick 2086 if (DataTypes::getRank(tempShape)==0) {
240 matt 1319
241    
242 jfenwick 1796 // get the space for the data vector
243     int len1 = DataTypes::noValues(tempShape);
244     DataVector temp_data(len1, 0.0, len1);
245 jfenwick 2105 temp_data.copyFromNumArray(asNumArray,1);
246 jfenwick 1796
247     int len = DataTypes::noValues(other.getDataPointShape());
248    
249     DataVector temp2_data(len, temp_data[0]/*temp_dataView()*/, len);
250     DataConstant* t=new DataConstant(other.getFunctionSpace(),other.getDataPointShape(),temp2_data);
251 jfenwick 1872 m_data=DataAbstract_ptr(t);
252 jfenwick 1796
253 jgs 94 } else {
254     //
255     // Create a DataConstant with the same sample shape as other
256 jfenwick 1796 DataConstant* t=new DataConstant(asNumArray,other.getFunctionSpace());
257 jfenwick 1872 // boost::shared_ptr<DataAbstract> sp(t);
258     // m_data=sp;
259     m_data=DataAbstract_ptr(t);
260 jgs 94 }
261 gross 783 m_protected=false;
262 jgs 94 }
263    
264 jgs 151 Data::~Data()
265     {
266    
267     }
268    
269 jfenwick 1796
270    
271     void
272     Data::initialise(const boost::python::numeric::array& value,
273     const FunctionSpace& what,
274     bool expanded)
275     {
276     //
277     // Construct a Data object of the appropriate type.
278     // Construct the object first as there seems to be a bug which causes
279     // undefined behaviour if an exception is thrown during construction
280     // within the shared_ptr constructor.
281     if (expanded) {
282     DataAbstract* temp=new DataExpanded(value, what);
283 jfenwick 1872 // boost::shared_ptr<DataAbstract> temp_data(temp);
284     // m_data=temp_data;
285     m_data=temp->getPtr();
286 jfenwick 1796 } else {
287     DataAbstract* temp=new DataConstant(value, what);
288 jfenwick 1872 // boost::shared_ptr<DataAbstract> temp_data(temp);
289     // m_data=temp_data;
290     m_data=temp->getPtr();
291 jfenwick 1796 }
292     }
293    
294    
295     void
296     Data::initialise(const DataTypes::ValueType& value,
297     const DataTypes::ShapeType& shape,
298     const FunctionSpace& what,
299     bool expanded)
300     {
301     //
302     // Construct a Data object of the appropriate type.
303     // Construct the object first as there seems to be a bug which causes
304     // undefined behaviour if an exception is thrown during construction
305     // within the shared_ptr constructor.
306     if (expanded) {
307     DataAbstract* temp=new DataExpanded(what, shape, value);
308 jfenwick 1872 // boost::shared_ptr<DataAbstract> temp_data(temp);
309     // m_data=temp_data;
310     m_data=temp->getPtr();
311 jfenwick 1796 } else {
312     DataAbstract* temp=new DataConstant(what, shape, value);
313 jfenwick 1872 // boost::shared_ptr<DataAbstract> temp_data(temp);
314     // m_data=temp_data;
315     m_data=temp->getPtr();
316 jfenwick 1796 }
317     }
318    
319    
320     // void
321     // Data::CompareDebug(const Data& rd)
322     // {
323     // using namespace std;
324     // bool mismatch=false;
325     // std::cout << "Comparing left and right" << endl;
326     // const DataTagged* left=dynamic_cast<DataTagged*>(m_data.get());
327     // const DataTagged* right=dynamic_cast<DataTagged*>(rd.m_data.get());
328     //
329     // if (left==0)
330     // {
331     // cout << "left arg is not a DataTagged\n";
332     // return;
333     // }
334     //
335     // if (right==0)
336     // {
337     // cout << "right arg is not a DataTagged\n";
338     // return;
339     // }
340     // cout << "Num elements=" << left->getVector().size() << ":" << right->getVector().size() << std::endl;
341     // cout << "Shapes ";
342     // if (left->getShape()==right->getShape())
343     // {
344     // cout << "ok\n";
345     // }
346     // else
347     // {
348     // cout << "Problem: shapes do not match\n";
349     // mismatch=true;
350     // }
351     // int lim=left->getVector().size();
352     // if (right->getVector().size()) lim=right->getVector().size();
353     // for (int i=0;i<lim;++i)
354     // {
355     // if (left->getVector()[i]!=right->getVector()[i])
356     // {
357     // cout << "[" << i << "] value mismatch " << left->getVector()[i] << ":" << right->getVector()[i] << endl;
358     // mismatch=true;
359     // }
360     // }
361     //
362     // // still need to check the tag map
363     // // also need to watch what is happening to function spaces, are they copied or what?
364     //
365     // const DataTagged::DataMapType& mapleft=left->getTagLookup();
366     // const DataTagged::DataMapType& mapright=right->getTagLookup();
367     //
368     // if (mapleft.size()!=mapright.size())
369     // {
370     // cout << "Maps are different sizes " << mapleft.size() << ":" << mapright.size() << endl;
371     // mismatch=true;
372     // cout << "Left map\n";
373     // DataTagged::DataMapType::const_iterator i,j;
374     // for (i=mapleft.begin();i!=mapleft.end();++i) {
375     // cout << "(" << i->first << "=>" << i->second << ")\n";
376     // }
377     // cout << "Right map\n";
378     // for (i=mapright.begin();i!=mapright.end();++i) {
379     // cout << "(" << i->first << "=>" << i->second << ")\n";
380     // }
381     // cout << "End map\n";
382     //
383     // }
384     //
385     // DataTagged::DataMapType::const_iterator i,j;
386     // for (i=mapleft.begin(),j=mapright.begin();i!=mapleft.end() && j!=mapright.end();++i,++j) {
387     // if ((i->first!=j->first) || (i->second!=j->second))
388     // {
389     // cout << "(" << i->first << "=>" << i->second << ")";
390     // cout << ":(" << j->first << "=>" << j->second << ") ";
391     // mismatch=true;
392     // }
393     // }
394     // if (mismatch)
395     // {
396     // cout << "#Mismatch\n";
397     // }
398     // }
399    
400 jgs 94 escriptDataC
401     Data::getDataC()
402     {
403     escriptDataC temp;
404     temp.m_dataPtr=(void*)this;
405     return temp;
406     }
407    
408     escriptDataC
409     Data::getDataC() const
410     {
411     escriptDataC temp;
412     temp.m_dataPtr=(void*)this;
413     return temp;
414     }
415    
416 jgs 121 const boost::python::tuple
417 jgs 94 Data::getShapeTuple() const
418     {
419 jfenwick 1796 const DataTypes::ShapeType& shape=getDataPointShape();
420 jgs 94 switch(getDataPointRank()) {
421     case 0:
422     return make_tuple();
423     case 1:
424     return make_tuple(long_(shape[0]));
425     case 2:
426     return make_tuple(long_(shape[0]),long_(shape[1]));
427     case 3:
428     return make_tuple(long_(shape[0]),long_(shape[1]),long_(shape[2]));
429     case 4:
430     return make_tuple(long_(shape[0]),long_(shape[1]),long_(shape[2]),long_(shape[3]));
431     default:
432     throw DataException("Error - illegal Data rank.");
433     }
434     }
435 jfenwick 1799
436    
437     // The different name is needed because boost has trouble with overloaded functions.
438     // It can't work out what type the function is based soley on its name.
439     // There are ways to fix this involving creating function pointer variables for each form
440     // but there doesn't seem to be a need given that the methods have the same name from the python point of view
441 jfenwick 2105 Data
442 jfenwick 1799 Data::copySelf()
443     {
444     DataAbstract* temp=m_data->deepCopy();
445 jfenwick 2105 return Data(temp);
446 jfenwick 1799 }
447    
448 jgs 94 void
449     Data::copy(const Data& other)
450     {
451 jfenwick 1799 DataAbstract* temp=other.m_data->deepCopy();
452 jfenwick 1872 DataAbstract_ptr p=temp->getPtr();
453     m_data=p;
454 jgs 94 }
455    
456 gross 1118
457 jfenwick 2005 Data
458     Data::delay()
459     {
460     DataLazy* dl=new DataLazy(m_data);
461     return Data(dl);
462     }
463    
464 jgs 94 void
465 jfenwick 2005 Data::delaySelf()
466     {
467     if (!isLazy())
468     {
469     m_data=(new DataLazy(m_data))->getPtr();
470     }
471     }
472    
473     void
474 gross 1118 Data::setToZero()
475     {
476 jfenwick 1803 if (isEmpty())
477 gross 1118 {
478 jfenwick 1803 throw DataException("Error - Operations not permitted on instances of DataEmpty.");
479     }
480 jfenwick 2005 m_data->setToZero();
481 gross 1118 }
482    
483     void
484 jgs 94 Data::copyWithMask(const Data& other,
485     const Data& mask)
486     {
487 jfenwick 1855 // 1. Interpolate if required so all Datas use the same FS as this
488     // 2. Tag or Expand so that all Data's are the same type
489     // 3. Iterate over the data vectors copying values where mask is >0
490 jfenwick 1803 if (other.isEmpty() || mask.isEmpty())
491     {
492     throw DataException("Error - copyWithMask not permitted using instances of DataEmpty.");
493     }
494 jfenwick 1855 Data other2(other);
495     Data mask2(mask);
496 jfenwick 2005 other2.resolve();
497     mask2.resolve();
498     this->resolve();
499 jfenwick 1855 FunctionSpace myFS=getFunctionSpace();
500     FunctionSpace oFS=other2.getFunctionSpace();
501     FunctionSpace mFS=mask2.getFunctionSpace();
502     if (oFS!=myFS)
503     {
504     if (other2.probeInterpolation(myFS))
505     {
506     other2=other2.interpolate(myFS);
507     }
508     else
509     {
510     throw DataException("Error - copyWithMask: other FunctionSpace is not compatible with this one.");
511     }
512     }
513     if (mFS!=myFS)
514     {
515     if (mask2.probeInterpolation(myFS))
516     {
517     mask2=mask2.interpolate(myFS);
518     }
519     else
520     {
521     throw DataException("Error - copyWithMask: mask FunctionSpace is not compatible with this one.");
522     }
523     }
524     // Ensure that all args have the same type
525     if (this->isExpanded() || mask2.isExpanded() || other2.isExpanded())
526     {
527     this->expand();
528     other2.expand();
529     mask2.expand();
530     }
531     else if (this->isTagged() || mask2.isTagged() || other2.isTagged())
532     {
533     this->tag();
534     other2.tag();
535     mask2.tag();
536     }
537     else if (this->isConstant() && mask2.isConstant() && other2.isConstant())
538     {
539     }
540     else
541     {
542     throw DataException("Error - Unknown DataAbstract passed to copyWithMask.");
543     }
544 jfenwick 2246 unsigned int selfrank=getDataPointRank();
545     unsigned int otherrank=other2.getDataPointRank();
546     unsigned int maskrank=mask2.getDataPointRank();
547     if ((selfrank==0) && (otherrank>0 || maskrank>0))
548     {
549     // to get here we must be copying from a large object into a scalar
550     // I am not allowing this.
551     // If you are calling copyWithMask then you are considering keeping some existing values
552     // and so I'm going to assume that you don't want your data objects getting a new shape.
553     throw DataException("Attempt to copyWithMask into a scalar from an object or mask with rank>0.");
554     }
555 jfenwick 1855 // Now we iterate over the elements
556 jfenwick 2005 DataVector& self=getReadyPtr()->getVector();
557     const DataVector& ovec=other2.getReadyPtr()->getVector();
558     const DataVector& mvec=mask2.getReadyPtr()->getVector();
559 jfenwick 2246
560     if ((selfrank>0) && (otherrank==0) &&(maskrank==0))
561 jfenwick 1855 {
562 jfenwick 2246 // Not allowing this combination.
563     // it is not clear what the rank of the target should be.
564     // Should it be filled with the scalar (rank stays the same);
565     // or should the target object be reshaped to be a scalar as well.
566     throw DataException("Attempt to copyWithMask from scalar mask and data into non-scalar target.");
567 jfenwick 1855 }
568 jfenwick 2246 if ((selfrank>0) && (otherrank>0) &&(maskrank==0))
569     {
570     if (mvec[0]>0) // copy whole object if scalar is >0
571     {
572     copy(other);
573     }
574     return;
575     }
576     if (isTagged()) // so all objects involved will also be tagged
577     {
578     // note the !
579     if (!((getDataPointShape()==mask2.getDataPointShape()) &&
580     ((other2.getDataPointShape()==mask2.getDataPointShape()) || (otherrank==0))))
581     {
582     throw DataException("copyWithMask, shape mismatch.");
583     }
584    
585     // We need to consider the possibility that tags are missing or in the wrong order
586     // My guiding assumption here is: All tagged Datas are assumed to have the default value for
587     // all tags which are not explicitly defined
588    
589     const DataTagged* mptr=dynamic_cast<const DataTagged*>(mask2.m_data.get());
590     const DataTagged* optr=dynamic_cast<const DataTagged*>(other2.m_data.get());
591     DataTagged* tptr=dynamic_cast<DataTagged*>(m_data.get());
592    
593     // first, add any tags missing from other or mask
594     const DataTagged::DataMapType& olookup=optr->getTagLookup();
595     const DataTagged::DataMapType& mlookup=mptr->getTagLookup();
596     DataTagged::DataMapType::const_iterator i; // i->first is a tag, i->second is an offset into memory
597     for (i=olookup.begin();i!=olookup.end();i++)
598     {
599     tptr->addTag(i->first);
600     }
601     for (i=mlookup.begin();i!=mlookup.end();i++) {
602     tptr->addTag(i->first);
603     }
604     // now we know that *this has all the required tags but they aren't guaranteed to be in
605     // the same order
606    
607     // There are two possibilities: 1. all objects have the same rank. 2. other is a scalar
608     if ((selfrank==otherrank) && (otherrank==maskrank))
609     {
610     for (i=olookup.begin();i!=olookup.end();i++)
611     {
612     // get the target offset
613     DataTypes::ValueType::size_type toff=tptr->getOffsetForTag(i->first);
614     DataTypes::ValueType::size_type moff=mptr->getOffsetForTag(i->first);
615     DataTypes::ValueType::size_type ooff=optr->getOffsetForTag(i->first);
616     for (int i=0;i<getDataPointSize();++i)
617     {
618     if (mvec[i+moff]>0)
619     {
620     self[i+toff]=ovec[i+ooff];
621     }
622     }
623     }
624     }
625     else // other is a scalar
626     {
627     for (i=mlookup.begin();i!=mlookup.end();i++)
628     {
629     // get the target offset
630     DataTypes::ValueType::size_type toff=tptr->getOffsetForTag(i->first);
631     DataTypes::ValueType::size_type moff=mptr->getOffsetForTag(i->first);
632     DataTypes::ValueType::size_type ooff=optr->getOffsetForTag(i->first);
633     for (int i=0;i<getDataPointSize();++i)
634     {
635     if (mvec[i+moff]>0)
636     {
637     self[i+toff]=ovec[ooff];
638     }
639     }
640     }
641     }
642    
643     return; // ugly
644     }
645     // mixed scalar and non-scalar operation
646     if ((selfrank>0) && (otherrank==0) && (mask2.getDataPointShape()==getDataPointShape()))
647     {
648     size_t num_points=self.size();
649     // OPENMP 3.0 allows unsigned loop vars.
650     #if defined(_OPENMP) && (_OPENMP < 200805)
651     long i;
652     #else
653     size_t i;
654     #endif
655     #pragma omp parallel for private(i) schedule(static)
656     for (i=0;i<num_points;++i)
657     {
658     if (mvec[i]>0)
659     {
660     self[i]=ovec[0];
661     }
662     }
663     return; // ugly!
664     }
665     // tagged data is already taken care of so we only need to worry about shapes
666     // special cases with scalars are already dealt with so all we need to worry about is shape
667     if ((getDataPointShape()!=other2.getDataPointShape()) || getDataPointShape()!=mask2.getDataPointShape())
668     {
669     ostringstream oss;
670     oss <<"Error - size mismatch in arguments to copyWithMask.";
671     oss << "\nself_shape=" << DataTypes::shapeToString(getDataPointShape());
672     oss << " other2_shape=" << DataTypes::shapeToString(other2.getDataPointShape());
673     oss << " mask2_shape=" << DataTypes::shapeToString(mask2.getDataPointShape());
674     throw DataException(oss.str());
675     }
676 jfenwick 1855 size_t num_points=self.size();
677 phornby 1921
678     // OPENMP 3.0 allows unsigned loop vars.
679     #if defined(_OPENMP) && (_OPENMP < 200805)
680     long i;
681     #else
682 phornby 1918 size_t i;
683 phornby 1921 #endif
684 jfenwick 1855 #pragma omp parallel for private(i) schedule(static)
685     for (i=0;i<num_points;++i)
686     {
687     if (mvec[i]>0)
688     {
689     self[i]=ovec[i];
690     }
691     }
692     }
693 jgs 94
694    
695    
696     bool
697     Data::isExpanded() const
698     {
699     DataExpanded* temp=dynamic_cast<DataExpanded*>(m_data.get());
700     return (temp!=0);
701     }
702    
703     bool
704     Data::isTagged() const
705     {
706     DataTagged* temp=dynamic_cast<DataTagged*>(m_data.get());
707     return (temp!=0);
708     }
709    
710     bool
711     Data::isEmpty() const
712     {
713     DataEmpty* temp=dynamic_cast<DataEmpty*>(m_data.get());
714     return (temp!=0);
715     }
716    
717     bool
718     Data::isConstant() const
719     {
720     DataConstant* temp=dynamic_cast<DataConstant*>(m_data.get());
721     return (temp!=0);
722     }
723    
724 jfenwick 2005 bool
725     Data::isLazy() const
726     {
727     return m_data->isLazy();
728     }
729    
730     // at the moment this is synonymous with !isLazy() but that could change
731     bool
732     Data::isReady() const
733     {
734     return (dynamic_cast<DataReady*>(m_data.get())!=0);
735     }
736    
737    
738 jgs 94 void
739 ksteube 1312 Data::setProtection()
740     {
741 gross 783 m_protected=true;
742     }
743    
744     bool
745 ksteube 1312 Data::isProtected() const
746     {
747 gross 783 return m_protected;
748     }
749    
750    
751    
752     void
753 jgs 94 Data::expand()
754     {
755     if (isConstant()) {
756     DataConstant* tempDataConst=dynamic_cast<DataConstant*>(m_data.get());
757     DataAbstract* temp=new DataExpanded(*tempDataConst);
758 jfenwick 1872 // shared_ptr<DataAbstract> temp_data(temp);
759     // m_data=temp_data;
760     m_data=temp->getPtr();
761 jgs 94 } else if (isTagged()) {
762     DataTagged* tempDataTag=dynamic_cast<DataTagged*>(m_data.get());
763     DataAbstract* temp=new DataExpanded(*tempDataTag);
764 jfenwick 1872 // shared_ptr<DataAbstract> temp_data(temp);
765     // m_data=temp_data;
766     m_data=temp->getPtr();
767 jgs 94 } else if (isExpanded()) {
768     //
769     // do nothing
770     } else if (isEmpty()) {
771     throw DataException("Error - Expansion of DataEmpty not possible.");
772 jfenwick 2005 } else if (isLazy()) {
773     resolve();
774     expand(); // resolve might not give us expanded data
775 jgs 94 } else {
776     throw DataException("Error - Expansion not implemented for this Data type.");
777     }
778     }
779    
780     void
781     Data::tag()
782     {
783     if (isConstant()) {
784     DataConstant* tempDataConst=dynamic_cast<DataConstant*>(m_data.get());
785     DataAbstract* temp=new DataTagged(*tempDataConst);
786 jfenwick 1872 // shared_ptr<DataAbstract> temp_data(temp);
787     // m_data=temp_data;
788     m_data=temp->getPtr();
789 jgs 94 } else if (isTagged()) {
790     // do nothing
791     } else if (isExpanded()) {
792     throw DataException("Error - Creating tag data from DataExpanded not possible.");
793     } else if (isEmpty()) {
794     throw DataException("Error - Creating tag data from DataEmpty not possible.");
795 jfenwick 2005 } else if (isLazy()) {
796     DataAbstract_ptr res=m_data->resolve();
797     if (m_data->isExpanded())
798     {
799     throw DataException("Error - data would resolve to DataExpanded, tagging is not possible.");
800     }
801     m_data=res;
802     tag();
803 jgs 94 } else {
804     throw DataException("Error - Tagging not implemented for this Data type.");
805     }
806     }
807    
808 jfenwick 2005 void
809     Data::resolve()
810     {
811     if (isLazy())
812     {
813     m_data=m_data->resolve();
814     }
815     }
816    
817    
818 gross 854 Data
819     Data::oneOver() const
820 jgs 102 {
821 jfenwick 2146 MAKELAZYOP(RECIP)
822 matt 1334 return C_TensorUnaryOperation(*this, bind1st(divides<double>(),1.));
823 jgs 102 }
824    
825 jgs 94 Data
826 gross 698 Data::wherePositive() const
827 jgs 94 {
828 jfenwick 2146 MAKELAZYOP(GZ)
829 matt 1334 return C_TensorUnaryOperation(*this, bind2nd(greater<double>(),0.0));
830 jgs 94 }
831    
832     Data
833 gross 698 Data::whereNegative() const
834 jgs 102 {
835 jfenwick 2146 MAKELAZYOP(LZ)
836 matt 1334 return C_TensorUnaryOperation(*this, bind2nd(less<double>(),0.0));
837 jgs 102 }
838    
839     Data
840 gross 698 Data::whereNonNegative() const
841 jgs 94 {
842 jfenwick 2146 MAKELAZYOP(GEZ)
843 matt 1334 return C_TensorUnaryOperation(*this, bind2nd(greater_equal<double>(),0.0));
844 jgs 94 }
845    
846     Data
847 gross 698 Data::whereNonPositive() const
848 jgs 94 {
849 jfenwick 2146 MAKELAZYOP(LEZ)
850 matt 1334 return C_TensorUnaryOperation(*this, bind2nd(less_equal<double>(),0.0));
851 jgs 94 }
852    
853     Data
854 jgs 571 Data::whereZero(double tol) const
855 jgs 94 {
856 jfenwick 2147 // Data dataAbs=abs();
857     // return C_TensorUnaryOperation(dataAbs, bind2nd(less_equal<double>(),tol));
858     MAKELAZYOPOFF(EZ,tol)
859     return C_TensorUnaryOperation(*this, bind2nd(AbsLTE(),tol));
860    
861 jgs 94 }
862    
863     Data
864 jgs 571 Data::whereNonZero(double tol) const
865 jgs 102 {
866 jfenwick 2147 // Data dataAbs=abs();
867     // return C_TensorUnaryOperation(dataAbs, bind2nd(greater<double>(),tol));
868     MAKELAZYOPOFF(NEZ,tol)
869     return C_TensorUnaryOperation(*this, bind2nd(AbsGT(),tol));
870    
871 jgs 102 }
872    
873     Data
874 jgs 94 Data::interpolate(const FunctionSpace& functionspace) const
875     {
876     return Data(*this,functionspace);
877     }
878    
879     bool
880     Data::probeInterpolation(const FunctionSpace& functionspace) const
881     {
882 jfenwick 2005 return getFunctionSpace().probeInterpolation(functionspace);
883 jgs 94 }
884    
885     Data
886     Data::gradOn(const FunctionSpace& functionspace) const
887     {
888 jfenwick 1803 if (isEmpty())
889     {
890     throw DataException("Error - operation not permitted on instances of DataEmpty.");
891     }
892 matt 1353 double blocktimer_start = blocktimer_time();
893 jgs 94 if (functionspace.getDomain()!=getDomain())
894     throw DataException("Error - gradient cannot be calculated on different domains.");
895 jfenwick 1796 DataTypes::ShapeType grad_shape=getDataPointShape();
896 jgs 94 grad_shape.push_back(functionspace.getDim());
897     Data out(0.0,grad_shape,functionspace,true);
898 jfenwick 1872 getDomain()->setToGradient(out,*this);
899 matt 1353 blocktimer_increment("grad()", blocktimer_start);
900 jgs 94 return out;
901     }
902    
903     Data
904     Data::grad() const
905     {
906 jfenwick 1803 if (isEmpty())
907     {
908     throw DataException("Error - operation not permitted on instances of DataEmpty.");
909     }
910 jfenwick 1872 return gradOn(escript::function(*getDomain()));
911 jgs 94 }
912    
913     int
914     Data::getDataPointSize() const
915     {
916 jfenwick 1796 return m_data->getNoValues();
917 jgs 94 }
918    
919 jfenwick 1796 DataTypes::ValueType::size_type
920 jgs 94 Data::getLength() const
921     {
922     return m_data->getLength();
923     }
924    
925 ksteube 1312 const
926 jgs 121 boost::python::numeric::array
927 ksteube 1312 Data:: getValueOfDataPoint(int dataPointNo)
928 jgs 121 {
929 gross 921 int i, j, k, l;
930 jfenwick 2005
931     FORCERESOLVE;
932    
933 jgs 121 //
934     // determine the rank and shape of each data point
935     int dataPointRank = getDataPointRank();
936 jfenwick 1796 const DataTypes::ShapeType& dataPointShape = getDataPointShape();
937 jgs 121
938     //
939     // create the numeric array to be returned
940     boost::python::numeric::array numArray(0.0);
941    
942     //
943 gross 921 // the shape of the returned numeric array will be the same
944     // as that of the data point
945     int arrayRank = dataPointRank;
946 jfenwick 1796 const DataTypes::ShapeType& arrayShape = dataPointShape;
947 jgs 121
948     //
949     // resize the numeric array to the shape just calculated
950 gross 921 if (arrayRank==0) {
951     numArray.resize(1);
952     }
953 jgs 121 if (arrayRank==1) {
954     numArray.resize(arrayShape[0]);
955     }
956     if (arrayRank==2) {
957     numArray.resize(arrayShape[0],arrayShape[1]);
958     }
959     if (arrayRank==3) {
960     numArray.resize(arrayShape[0],arrayShape[1],arrayShape[2]);
961     }
962     if (arrayRank==4) {
963     numArray.resize(arrayShape[0],arrayShape[1],arrayShape[2],arrayShape[3]);
964     }
965    
966 gross 921 if (getNumDataPointsPerSample()>0) {
967     int sampleNo = dataPointNo/getNumDataPointsPerSample();
968     int dataPointNoInSample = dataPointNo - sampleNo * getNumDataPointsPerSample();
969     //
970     // Check a valid sample number has been supplied
971 trankine 924 if ((sampleNo >= getNumSamples()) || (sampleNo < 0 )) {
972 gross 921 throw DataException("Error - Data::convertToNumArray: invalid sampleNo.");
973     }
974 ksteube 1312
975 gross 921 //
976     // Check a valid data point number has been supplied
977 trankine 924 if ((dataPointNoInSample >= getNumDataPointsPerSample()) || (dataPointNoInSample < 0)) {
978 gross 921 throw DataException("Error - Data::convertToNumArray: invalid dataPointNoInSample.");
979     }
980     // TODO: global error handling
981     // create a view of the data if it is stored locally
982 jfenwick 1796 // DataArrayView dataPointView = getDataPoint(sampleNo, dataPointNoInSample);
983     DataTypes::ValueType::size_type offset=getDataOffset(sampleNo, dataPointNoInSample);
984 ksteube 1312
985 jfenwick 1796
986 gross 921 switch( dataPointRank ){
987     case 0 :
988 jfenwick 1796 numArray[0] = getDataAtOffset(offset);
989 gross 921 break;
990 ksteube 1312 case 1 :
991 gross 921 for( i=0; i<dataPointShape[0]; i++ )
992 jfenwick 1796 numArray[i]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i));
993 gross 921 break;
994 ksteube 1312 case 2 :
995 gross 921 for( i=0; i<dataPointShape[0]; i++ )
996     for( j=0; j<dataPointShape[1]; j++)
997 jfenwick 1796 numArray[make_tuple(i,j)]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i,j));
998 gross 921 break;
999 ksteube 1312 case 3 :
1000 gross 921 for( i=0; i<dataPointShape[0]; i++ )
1001     for( j=0; j<dataPointShape[1]; j++ )
1002     for( k=0; k<dataPointShape[2]; k++)
1003 jfenwick 1796 numArray[make_tuple(i,j,k)]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i,j,k));
1004 gross 921 break;
1005     case 4 :
1006     for( i=0; i<dataPointShape[0]; i++ )
1007     for( j=0; j<dataPointShape[1]; j++ )
1008     for( k=0; k<dataPointShape[2]; k++ )
1009     for( l=0; l<dataPointShape[3]; l++)
1010 jfenwick 1796 numArray[make_tuple(i,j,k,l)]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i,j,k,l));
1011 gross 921 break;
1012     }
1013 jgs 117 }
1014     //
1015 gross 921 // return the array
1016 jgs 117 return numArray;
1017     }
1018 jfenwick 1796
1019 gross 921 void
1020 ksteube 1312 Data::setValueOfDataPointToPyObject(int dataPointNo, const boost::python::object& py_object)
1021 jgs 121 {
1022 gross 1034 // this will throw if the value cannot be represented
1023     boost::python::numeric::array num_array(py_object);
1024     setValueOfDataPointToArray(dataPointNo,num_array);
1025     }
1026    
1027     void
1028     Data::setValueOfDataPointToArray(int dataPointNo, const boost::python::numeric::array& num_array)
1029     {
1030 gross 921 if (isProtected()) {
1031     throw DataException("Error - attempt to update protected Data object.");
1032     }
1033 jfenwick 2005 FORCERESOLVE;
1034 gross 921 //
1035     // check rank
1036 jfenwick 1977 if (static_cast<unsigned int>(num_array.getrank())<getDataPointRank())
1037 gross 921 throw DataException("Rank of numarray does not match Data object rank");
1038 bcumming 790
1039 jgs 121 //
1040 gross 921 // check shape of num_array
1041 phornby 1953 for (unsigned int i=0; i<getDataPointRank(); i++) {
1042 gross 921 if (extract<int>(num_array.getshape()[i])!=getDataPointShape()[i])
1043     throw DataException("Shape of numarray does not match Data object rank");
1044 jgs 121 }
1045     //
1046 gross 921 // make sure data is expanded:
1047 gross 1859 //
1048 gross 921 if (!isExpanded()) {
1049     expand();
1050 jgs 121 }
1051 gross 921 if (getNumDataPointsPerSample()>0) {
1052     int sampleNo = dataPointNo/getNumDataPointsPerSample();
1053     int dataPointNoInSample = dataPointNo - sampleNo * getNumDataPointsPerSample();
1054     m_data->copyToDataPoint(sampleNo, dataPointNoInSample,num_array);
1055     } else {
1056     m_data->copyToDataPoint(-1, 0,num_array);
1057     }
1058     }
1059 jgs 121
1060 gross 922 void
1061     Data::setValueOfDataPoint(int dataPointNo, const double value)
1062     {
1063     if (isProtected()) {
1064     throw DataException("Error - attempt to update protected Data object.");
1065     }
1066     //
1067     // make sure data is expanded:
1068 jfenwick 2005 FORCERESOLVE;
1069 gross 922 if (!isExpanded()) {
1070     expand();
1071     }
1072     if (getNumDataPointsPerSample()>0) {
1073     int sampleNo = dataPointNo/getNumDataPointsPerSample();
1074     int dataPointNoInSample = dataPointNo - sampleNo * getNumDataPointsPerSample();
1075     m_data->copyToDataPoint(sampleNo, dataPointNoInSample,value);
1076     } else {
1077     m_data->copyToDataPoint(-1, 0,value);
1078     }
1079     }
1080    
1081 ksteube 1312 const
1082 gross 921 boost::python::numeric::array
1083 ksteube 1312 Data::getValueOfGlobalDataPoint(int procNo, int dataPointNo)
1084 gross 921 {
1085     size_t length=0;
1086     int i, j, k, l, pos;
1087 jfenwick 2005 FORCERESOLVE;
1088 jgs 121 //
1089     // determine the rank and shape of each data point
1090     int dataPointRank = getDataPointRank();
1091 jfenwick 1796 const DataTypes::ShapeType& dataPointShape = getDataPointShape();
1092 jgs 121
1093     //
1094     // create the numeric array to be returned
1095     boost::python::numeric::array numArray(0.0);
1096    
1097     //
1098     // the shape of the returned numeric array will be the same
1099     // as that of the data point
1100     int arrayRank = dataPointRank;
1101 jfenwick 1796 const DataTypes::ShapeType& arrayShape = dataPointShape;
1102 jgs 121
1103     //
1104     // resize the numeric array to the shape just calculated
1105     if (arrayRank==0) {
1106     numArray.resize(1);
1107     }
1108     if (arrayRank==1) {
1109     numArray.resize(arrayShape[0]);
1110     }
1111     if (arrayRank==2) {
1112     numArray.resize(arrayShape[0],arrayShape[1]);
1113     }
1114     if (arrayRank==3) {
1115     numArray.resize(arrayShape[0],arrayShape[1],arrayShape[2]);
1116     }
1117     if (arrayRank==4) {
1118     numArray.resize(arrayShape[0],arrayShape[1],arrayShape[2],arrayShape[3]);
1119     }
1120    
1121 gross 921 // added for the MPI communication
1122     length=1;
1123     for( i=0; i<arrayRank; i++ ) length *= arrayShape[i];
1124     double *tmpData = new double[length];
1125 bcumming 790
1126 jgs 121 //
1127     // load the values for the data point into the numeric array.
1128 bcumming 790
1129     // updated for the MPI case
1130     if( get_MPIRank()==procNo ){
1131 gross 921 if (getNumDataPointsPerSample()>0) {
1132     int sampleNo = dataPointNo/getNumDataPointsPerSample();
1133     int dataPointNoInSample = dataPointNo - sampleNo * getNumDataPointsPerSample();
1134     //
1135     // Check a valid sample number has been supplied
1136 trankine 924 if ((sampleNo >= getNumSamples()) || (sampleNo < 0 )) {
1137 gross 921 throw DataException("Error - Data::convertToNumArray: invalid sampleNo.");
1138     }
1139 ksteube 1312
1140 gross 921 //
1141     // Check a valid data point number has been supplied
1142 trankine 924 if ((dataPointNoInSample >= getNumDataPointsPerSample()) || (dataPointNoInSample < 0)) {
1143 gross 921 throw DataException("Error - Data::convertToNumArray: invalid dataPointNoInSample.");
1144     }
1145     // TODO: global error handling
1146 bcumming 790 // create a view of the data if it is stored locally
1147 jfenwick 1796 //DataArrayView dataPointView = getDataPoint(sampleNo, dataPointNoInSample);
1148     DataTypes::ValueType::size_type offset=getDataOffset(sampleNo, dataPointNoInSample);
1149 ksteube 1312
1150 bcumming 790 // pack the data from the view into tmpData for MPI communication
1151     pos=0;
1152     switch( dataPointRank ){
1153     case 0 :
1154 jfenwick 1796 tmpData[0] = getDataAtOffset(offset);
1155 bcumming 790 break;
1156 ksteube 1312 case 1 :
1157 bcumming 790 for( i=0; i<dataPointShape[0]; i++ )
1158 jfenwick 1796 tmpData[i]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i));
1159 bcumming 790 break;
1160 ksteube 1312 case 2 :
1161 bcumming 790 for( i=0; i<dataPointShape[0]; i++ )
1162     for( j=0; j<dataPointShape[1]; j++, pos++ )
1163 jfenwick 1796 tmpData[pos]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i,j));
1164 bcumming 790 break;
1165 ksteube 1312 case 3 :
1166 bcumming 790 for( i=0; i<dataPointShape[0]; i++ )
1167     for( j=0; j<dataPointShape[1]; j++ )
1168     for( k=0; k<dataPointShape[2]; k++, pos++ )
1169 jfenwick 1796 tmpData[pos]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i,j,k));
1170 bcumming 790 break;
1171     case 4 :
1172     for( i=0; i<dataPointShape[0]; i++ )
1173     for( j=0; j<dataPointShape[1]; j++ )
1174     for( k=0; k<dataPointShape[2]; k++ )
1175     for( l=0; l<dataPointShape[3]; l++, pos++ )
1176 jfenwick 1796 tmpData[pos]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i,j,k,l));
1177 bcumming 790 break;
1178     }
1179 gross 921 }
1180 bcumming 790 }
1181 ksteube 1312 #ifdef PASO_MPI
1182 gross 921 // broadcast the data to all other processes
1183     MPI_Bcast( tmpData, length, MPI_DOUBLE, procNo, get_MPIComm() );
1184     #endif
1185 bcumming 790
1186     // unpack the data
1187     switch( dataPointRank ){
1188     case 0 :
1189 gross 921 numArray[0]=tmpData[0];
1190 bcumming 790 break;
1191 ksteube 1312 case 1 :
1192 bcumming 790 for( i=0; i<dataPointShape[0]; i++ )
1193     numArray[i]=tmpData[i];
1194     break;
1195 ksteube 1312 case 2 :
1196 bcumming 790 for( i=0; i<dataPointShape[0]; i++ )
1197     for( j=0; j<dataPointShape[1]; j++ )
1198 gross 921 numArray[make_tuple(i,j)]=tmpData[i+j*dataPointShape[0]];
1199 bcumming 790 break;
1200 ksteube 1312 case 3 :
1201 bcumming 790 for( i=0; i<dataPointShape[0]; i++ )
1202     for( j=0; j<dataPointShape[1]; j++ )
1203     for( k=0; k<dataPointShape[2]; k++ )
1204 gross 921 numArray[make_tuple(i,j,k)]=tmpData[i+dataPointShape[0]*(j*+k*dataPointShape[1])];
1205 bcumming 790 break;
1206     case 4 :
1207     for( i=0; i<dataPointShape[0]; i++ )
1208     for( j=0; j<dataPointShape[1]; j++ )
1209     for( k=0; k<dataPointShape[2]; k++ )
1210     for( l=0; l<dataPointShape[3]; l++ )
1211 gross 921 numArray[make_tuple(i,j,k,l)]=tmpData[i+dataPointShape[0]*(j*+dataPointShape[1]*(k+l*dataPointShape[2]))];
1212 bcumming 790 break;
1213     }
1214    
1215 ksteube 1312 delete [] tmpData;
1216 jgs 121 //
1217     // return the loaded array
1218     return numArray;
1219     }
1220    
1221 gross 921
1222 jfenwick 2005 boost::python::numeric::array
1223     Data::integrate_const() const
1224     {
1225     if (isLazy())
1226     {
1227     throw DataException("Error - cannot integrate for constant lazy data.");
1228     }
1229     return integrateWorker();
1230     }
1231 gross 921
1232 jgs 121 boost::python::numeric::array
1233 jfenwick 2005 Data::integrate()
1234 jgs 94 {
1235 jfenwick 2005 if (isLazy())
1236     {
1237     expand();
1238     }
1239     return integrateWorker();
1240     }
1241    
1242    
1243    
1244     boost::python::numeric::array
1245     Data::integrateWorker() const
1246     {
1247 jgs 94 int index;
1248     int rank = getDataPointRank();
1249 jfenwick 1796 DataTypes::ShapeType shape = getDataPointShape();
1250 ksteube 1312 int dataPointSize = getDataPointSize();
1251 jgs 94
1252     //
1253     // calculate the integral values
1254 ksteube 1312 vector<double> integrals(dataPointSize);
1255     vector<double> integrals_local(dataPointSize);
1256 jfenwick 2089 const AbstractContinuousDomain* dom=dynamic_cast<const AbstractContinuousDomain*>(getDomain().get());
1257     if (dom==0)
1258     {
1259     throw DataException("Can not integrate over non-continuous domains.");
1260     }
1261 ksteube 1312 #ifdef PASO_MPI
1262 jfenwick 2089 dom->setToIntegrals(integrals_local,*this);
1263 ksteube 1312 // Global sum: use an array instead of a vector because elements of array are guaranteed to be contiguous in memory
1264     double *tmp = new double[dataPointSize];
1265     double *tmp_local = new double[dataPointSize];
1266     for (int i=0; i<dataPointSize; i++) { tmp_local[i] = integrals_local[i]; }
1267     MPI_Allreduce( &tmp_local[0], &tmp[0], dataPointSize, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD );
1268     for (int i=0; i<dataPointSize; i++) { integrals[i] = tmp[i]; }
1269     delete[] tmp;
1270     delete[] tmp_local;
1271     #else
1272 jfenwick 2089 dom->setToIntegrals(integrals,*this);
1273 ksteube 1312 #endif
1274 jgs 94
1275     //
1276     // create the numeric array to be returned
1277     // and load the array with the integral values
1278     boost::python::numeric::array bp_array(1.0);
1279     if (rank==0) {
1280 jgs 108 bp_array.resize(1);
1281 jgs 94 index = 0;
1282     bp_array[0] = integrals[index];
1283     }
1284     if (rank==1) {
1285     bp_array.resize(shape[0]);
1286     for (int i=0; i<shape[0]; i++) {
1287     index = i;
1288     bp_array[i] = integrals[index];
1289     }
1290     }
1291     if (rank==2) {
1292 gross 436 bp_array.resize(shape[0],shape[1]);
1293     for (int i=0; i<shape[0]; i++) {
1294     for (int j=0; j<shape[1]; j++) {
1295     index = i + shape[0] * j;
1296     bp_array[make_tuple(i,j)] = integrals[index];
1297     }
1298     }
1299 jgs 94 }
1300     if (rank==3) {
1301     bp_array.resize(shape[0],shape[1],shape[2]);
1302     for (int i=0; i<shape[0]; i++) {
1303     for (int j=0; j<shape[1]; j++) {
1304     for (int k=0; k<shape[2]; k++) {
1305     index = i + shape[0] * ( j + shape[1] * k );
1306 gross 436 bp_array[make_tuple(i,j,k)] = integrals[index];
1307 jgs 94 }
1308     }
1309     }
1310     }
1311     if (rank==4) {
1312     bp_array.resize(shape[0],shape[1],shape[2],shape[3]);
1313     for (int i=0; i<shape[0]; i++) {
1314     for (int j=0; j<shape[1]; j++) {
1315     for (int k=0; k<shape[2]; k++) {
1316     for (int l=0; l<shape[3]; l++) {
1317     index = i + shape[0] * ( j + shape[1] * ( k + shape[2] * l ) );
1318 gross 436 bp_array[make_tuple(i,j,k,l)] = integrals[index];
1319 jgs 94 }
1320     }
1321     }
1322     }
1323     }
1324    
1325     //
1326     // return the loaded array
1327     return bp_array;
1328     }
1329    
1330     Data
1331     Data::sin() const
1332     {
1333 jfenwick 2146 MAKELAZYOP(SIN)
1334 matt 1349 return C_TensorUnaryOperation<double (*)(double)>(*this, ::sin);
1335 jgs 94 }
1336    
1337     Data
1338     Data::cos() const
1339     {
1340 jfenwick 2146 MAKELAZYOP(COS)
1341 matt 1349 return C_TensorUnaryOperation<double (*)(double)>(*this, ::cos);
1342 jgs 94 }
1343    
1344     Data
1345     Data::tan() const
1346     {
1347 jfenwick 2146 MAKELAZYOP(TAN)
1348 matt 1349 return C_TensorUnaryOperation<double (*)(double)>(*this, ::tan);
1349 jgs 94 }
1350    
1351     Data
1352 jgs 150 Data::asin() const
1353     {
1354 jfenwick 2146 MAKELAZYOP(ASIN)
1355 matt 1350 return C_TensorUnaryOperation<double (*)(double)>(*this, ::asin);
1356 jgs 150 }
1357    
1358     Data
1359     Data::acos() const
1360     {
1361 jfenwick 2146 MAKELAZYOP(ACOS)
1362 matt 1350 return C_TensorUnaryOperation<double (*)(double)>(*this, ::acos);
1363 jgs 150 }
1364    
1365 phornby 1026
1366 jgs 150 Data
1367     Data::atan() const
1368     {
1369 jfenwick 2146 MAKELAZYOP(ATAN)
1370 matt 1350 return C_TensorUnaryOperation<double (*)(double)>(*this, ::atan);
1371 jgs 150 }
1372    
1373     Data
1374     Data::sinh() const
1375     {
1376 jfenwick 2146 MAKELAZYOP(SINH)
1377 jfenwick 2005 return C_TensorUnaryOperation<double (*)(double)>(*this, ::sinh);
1378 jgs 150 }
1379    
1380     Data
1381     Data::cosh() const
1382     {
1383 jfenwick 2146 MAKELAZYOP(COSH)
1384 jfenwick 2005 return C_TensorUnaryOperation<double (*)(double)>(*this, ::cosh);
1385 jgs 150 }
1386    
1387     Data
1388     Data::tanh() const
1389     {
1390 jfenwick 2146 MAKELAZYOP(TANH)
1391 jfenwick 2005 return C_TensorUnaryOperation<double (*)(double)>(*this, ::tanh);
1392 jgs 150 }
1393    
1394 phornby 1026
1395 jgs 150 Data
1396 phornby 1026 Data::erf() const
1397     {
1398 phornby 2048 #if defined (_WIN32) && !defined(__INTEL_COMPILER)
1399 gross 1028 throw DataException("Error - Data:: erf function is not supported on _WIN32 platforms.");
1400     #else
1401 jfenwick 2146 MAKELAZYOP(ERF)
1402 matt 1334 return C_TensorUnaryOperation(*this, ::erf);
1403 phornby 1026 #endif
1404     }
1405    
1406     Data
1407 jgs 150 Data::asinh() const
1408     {
1409 jfenwick 2146 MAKELAZYOP(ASINH)
1410 phornby 2048 #if defined (_WIN32) && !defined(__INTEL_COMPILER)
1411 matt 1334 return C_TensorUnaryOperation(*this, escript::asinh_substitute);
1412 phornby 1032 #else
1413 matt 1334 return C_TensorUnaryOperation(*this, ::asinh);
1414 phornby 1032 #endif
1415 jgs 150 }
1416    
1417     Data
1418     Data::acosh() const
1419     {
1420 jfenwick 2146 MAKELAZYOP(ACOSH)
1421 phornby 2048 #if defined (_WIN32) && !defined(__INTEL_COMPILER)
1422 matt 1334 return C_TensorUnaryOperation(*this, escript::acosh_substitute);
1423 phornby 1032 #else
1424 matt 1334 return C_TensorUnaryOperation(*this, ::acosh);
1425 phornby 1032 #endif
1426 jgs 150 }
1427    
1428     Data
1429     Data::atanh() const
1430     {
1431 jfenwick 2146 MAKELAZYOP(ATANH)
1432 phornby 2048 #if defined (_WIN32) && !defined(__INTEL_COMPILER)
1433 matt 1334 return C_TensorUnaryOperation(*this, escript::atanh_substitute);
1434 phornby 1032 #else
1435 matt 1334 return C_TensorUnaryOperation(*this, ::atanh);
1436 phornby 1032 #endif
1437 jgs 150 }
1438    
1439     Data
1440 gross 286 Data::log10() const
1441 jfenwick 2146 {
1442     MAKELAZYOP(LOG10)
1443 matt 1350 return C_TensorUnaryOperation<double (*)(double)>(*this, ::log10);
1444 jgs 94 }
1445    
1446     Data
1447 gross 286 Data::log() const
1448 jgs 94 {
1449 jfenwick 2146 MAKELAZYOP(LOG)
1450 matt 1350 return C_TensorUnaryOperation<double (*)(double)>(*this, ::log);
1451 jgs 94 }
1452    
1453 jgs 106 Data
1454     Data::sign() const
1455 jgs 94 {
1456 jfenwick 2146 MAKELAZYOP(SIGN)
1457 matt 1334 return C_TensorUnaryOperation(*this, escript::fsign);
1458 jgs 94 }
1459    
1460 jgs 106 Data
1461     Data::abs() const
1462 jgs 94 {
1463 jfenwick 2146 MAKELAZYOP(ABS)
1464 matt 1350 return C_TensorUnaryOperation<double (*)(double)>(*this, ::fabs);
1465 jgs 94 }
1466    
1467 jgs 106 Data
1468     Data::neg() const
1469 jgs 94 {
1470 jfenwick 2146 MAKELAZYOP(NEG)
1471 matt 1334 return C_TensorUnaryOperation(*this, negate<double>());
1472 jgs 94 }
1473    
1474 jgs 102 Data
1475 jgs 106 Data::pos() const
1476 jgs 94 {
1477 jfenwick 2005 // not doing lazy check here is deliberate.
1478     // since a deep copy of lazy data should be cheap, I'll just let it happen now
1479 jgs 148 Data result;
1480     // perform a deep copy
1481     result.copy(*this);
1482     return result;
1483 jgs 102 }
1484    
1485     Data
1486 jgs 106 Data::exp() const
1487 jfenwick 2146 {
1488     MAKELAZYOP(EXP)
1489 matt 1350 return C_TensorUnaryOperation<double (*)(double)>(*this, ::exp);
1490 jgs 102 }
1491    
1492     Data
1493 jgs 106 Data::sqrt() const
1494 jgs 102 {
1495 jfenwick 2146 MAKELAZYOP(SQRT)
1496 matt 1350 return C_TensorUnaryOperation<double (*)(double)>(*this, ::sqrt);
1497 jgs 102 }
1498    
1499 jgs 106 double
1500 jfenwick 2005 Data::Lsup_const() const
1501 jgs 102 {
1502 jfenwick 2005 if (isLazy())
1503     {
1504     throw DataException("Error - cannot compute Lsup for constant lazy data.");
1505     }
1506     return LsupWorker();
1507     }
1508    
1509     double
1510     Data::Lsup()
1511     {
1512     if (isLazy())
1513     {
1514 jfenwick 2085 resolve();
1515 jfenwick 2005 }
1516     return LsupWorker();
1517     }
1518    
1519     double
1520     Data::sup_const() const
1521     {
1522     if (isLazy())
1523     {
1524     throw DataException("Error - cannot compute sup for constant lazy data.");
1525     }
1526     return supWorker();
1527     }
1528    
1529     double
1530     Data::sup()
1531     {
1532     if (isLazy())
1533     {
1534 jfenwick 2085 resolve();
1535 jfenwick 2005 }
1536     return supWorker();
1537     }
1538    
1539     double
1540     Data::inf_const() const
1541     {
1542     if (isLazy())
1543     {
1544     throw DataException("Error - cannot compute inf for constant lazy data.");
1545     }
1546     return infWorker();
1547     }
1548    
1549     double
1550     Data::inf()
1551     {
1552     if (isLazy())
1553     {
1554 jfenwick 2085 resolve();
1555 jfenwick 2005 }
1556     return infWorker();
1557     }
1558    
1559     double
1560     Data::LsupWorker() const
1561     {
1562 phornby 1455 double localValue;
1563 jgs 106 //
1564     // set the initial absolute maximum value to zero
1565 bcumming 751
1566 jgs 147 AbsMax abs_max_func;
1567 bcumming 751 localValue = algorithm(abs_max_func,0);
1568     #ifdef PASO_MPI
1569 phornby 1455 double globalValue;
1570 bcumming 751 MPI_Allreduce( &localValue, &globalValue, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD );
1571     return globalValue;
1572     #else
1573     return localValue;
1574     #endif
1575 jgs 117 }
1576    
1577     double
1578 jfenwick 2005 Data::supWorker() const
1579 jgs 102 {
1580 phornby 1455 double localValue;
1581 jgs 106 //
1582     // set the initial maximum value to min possible double
1583 jgs 147 FMax fmax_func;
1584 bcumming 751 localValue = algorithm(fmax_func,numeric_limits<double>::max()*-1);
1585     #ifdef PASO_MPI
1586 phornby 1455 double globalValue;
1587 bcumming 751 MPI_Allreduce( &localValue, &globalValue, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD );
1588     return globalValue;
1589     #else
1590     return localValue;
1591     #endif
1592 jgs 102 }
1593    
1594 jgs 106 double
1595 jfenwick 2005 Data::infWorker() const
1596 jgs 102 {
1597 phornby 1455 double localValue;
1598 jgs 106 //
1599     // set the initial minimum value to max possible double
1600 jgs 147 FMin fmin_func;
1601 bcumming 751 localValue = algorithm(fmin_func,numeric_limits<double>::max());
1602     #ifdef PASO_MPI
1603 phornby 1455 double globalValue;
1604 bcumming 751 MPI_Allreduce( &localValue, &globalValue, 1, MPI_DOUBLE, MPI_MIN, MPI_COMM_WORLD );
1605     return globalValue;
1606     #else
1607     return localValue;
1608     #endif
1609 jgs 102 }
1610    
1611 bcumming 751 /* TODO */
1612     /* global reduction */
1613 jgs 102 Data
1614 jgs 106 Data::maxval() const
1615 jgs 102 {
1616 jfenwick 2037 if (isLazy())
1617     {
1618     Data temp(*this); // to get around the fact that you can't resolve a const Data
1619     temp.resolve();
1620     return temp.maxval();
1621     }
1622 jgs 113 //
1623     // set the initial maximum value to min possible double
1624 jgs 147 FMax fmax_func;
1625     return dp_algorithm(fmax_func,numeric_limits<double>::max()*-1);
1626 jgs 102 }
1627    
1628     Data
1629 jgs 106 Data::minval() const
1630 jgs 102 {
1631 jfenwick 2037 if (isLazy())
1632     {
1633     Data temp(*this); // to get around the fact that you can't resolve a const Data
1634     temp.resolve();
1635     return temp.minval();
1636     }
1637 jgs 113 //
1638     // set the initial minimum value to max possible double
1639 jgs 147 FMin fmin_func;
1640     return dp_algorithm(fmin_func,numeric_limits<double>::max());
1641 jgs 102 }
1642    
1643 jgs 123 Data
1644 gross 804 Data::swapaxes(const int axis0, const int axis1) const
1645 jgs 123 {
1646 jfenwick 2085 if (isLazy())
1647     {
1648     Data temp(*this);
1649     temp.resolve();
1650     return temp.swapaxes(axis0,axis1);
1651     }
1652 gross 804 int axis0_tmp,axis1_tmp;
1653 jfenwick 1796 DataTypes::ShapeType s=getDataPointShape();
1654     DataTypes::ShapeType ev_shape;
1655 gross 800 // Here's the equivalent of python s_out=s[axis_offset:]+s[:axis_offset]
1656     // which goes thru all shape vector elements starting with axis_offset (at index=rank wrap around to 0)
1657     int rank=getDataPointRank();
1658 gross 804 if (rank<2) {
1659     throw DataException("Error - Data::swapaxes argument must have at least rank 2.");
1660 gross 800 }
1661 gross 804 if (axis0<0 || axis0>rank-1) {
1662     throw DataException("Error - Data::swapaxes: axis0 must be between 0 and rank-1=" + rank-1);
1663     }
1664     if (axis1<0 || axis1>rank-1) {
1665     throw DataException("Error - Data::swapaxes: axis1 must be between 0 and rank-1=" + rank-1);
1666     }
1667     if (axis0 == axis1) {
1668     throw DataException("Error - Data::swapaxes: axis indices must be different.");
1669     }
1670     if (axis0 > axis1) {
1671     axis0_tmp=axis1;
1672     axis1_tmp=axis0;
1673     } else {
1674     axis0_tmp=axis0;
1675     axis1_tmp=axis1;
1676     }
1677 gross 800 for (int i=0; i<rank; i++) {
1678 gross 804 if (i == axis0_tmp) {
1679 ksteube 1312 ev_shape.push_back(s[axis1_tmp]);
1680 gross 804 } else if (i == axis1_tmp) {
1681 ksteube 1312 ev_shape.push_back(s[axis0_tmp]);
1682 gross 800 } else {
1683 ksteube 1312 ev_shape.push_back(s[i]);
1684 gross 800 }
1685     }
1686     Data ev(0.,ev_shape,getFunctionSpace());
1687     ev.typeMatchRight(*this);
1688 gross 804 m_data->swapaxes(ev.m_data.get(), axis0_tmp, axis1_tmp);
1689 gross 800 return ev;
1690    
1691 jgs 123 }
1692    
1693     Data
1694 ksteube 775 Data::symmetric() const
1695 jgs 123 {
1696 ksteube 775 // check input
1697 jfenwick 1796 DataTypes::ShapeType s=getDataPointShape();
1698 ksteube 775 if (getDataPointRank()==2) {
1699 ksteube 1312 if(s[0] != s[1])
1700 ksteube 775 throw DataException("Error - Data::symmetric can only be calculated for rank 2 object with equal first and second dimension.");
1701     }
1702     else if (getDataPointRank()==4) {
1703     if(!(s[0] == s[2] && s[1] == s[3]))
1704     throw DataException("Error - Data::symmetric can only be calculated for rank 4 object with dim0==dim2 and dim1==dim3.");
1705     }
1706     else {
1707     throw DataException("Error - Data::symmetric can only be calculated for rank 2 or 4 object.");
1708     }
1709 jfenwick 2146 MAKELAZYOP(SYM)
1710 ksteube 775 Data ev(0.,getDataPointShape(),getFunctionSpace());
1711     ev.typeMatchRight(*this);
1712     m_data->symmetric(ev.m_data.get());
1713     return ev;
1714     }
1715    
1716     Data
1717     Data::nonsymmetric() const
1718     {
1719 jfenwick 2146 MAKELAZYOP(NSYM)
1720 ksteube 775 // check input
1721 jfenwick 1796 DataTypes::ShapeType s=getDataPointShape();
1722 ksteube 775 if (getDataPointRank()==2) {
1723 ksteube 1312 if(s[0] != s[1])
1724 ksteube 775 throw DataException("Error - Data::nonsymmetric can only be calculated for rank 2 object with equal first and second dimension.");
1725 jfenwick 1796 DataTypes::ShapeType ev_shape;
1726 ksteube 775 ev_shape.push_back(s[0]);
1727     ev_shape.push_back(s[1]);
1728     Data ev(0.,ev_shape,getFunctionSpace());
1729     ev.typeMatchRight(*this);
1730     m_data->nonsymmetric(ev.m_data.get());
1731     return ev;
1732     }
1733     else if (getDataPointRank()==4) {
1734     if(!(s[0] == s[2] && s[1] == s[3]))
1735     throw DataException("Error - Data::nonsymmetric can only be calculated for rank 4 object with dim0==dim2 and dim1==dim3.");
1736 jfenwick 1796 DataTypes::ShapeType ev_shape;
1737 ksteube 775 ev_shape.push_back(s[0]);
1738     ev_shape.push_back(s[1]);
1739     ev_shape.push_back(s[2]);
1740     ev_shape.push_back(s[3]);
1741     Data ev(0.,ev_shape,getFunctionSpace());
1742     ev.typeMatchRight(*this);
1743     m_data->nonsymmetric(ev.m_data.get());
1744     return ev;
1745     }
1746     else {
1747     throw DataException("Error - Data::nonsymmetric can only be calculated for rank 2 or 4 object.");
1748     }
1749     }
1750    
1751     Data
1752 gross 800 Data::trace(int axis_offset) const
1753 jfenwick 2084 {
1754 jfenwick 2146 MAKELAZYOPOFF(TRACE,axis_offset)
1755 jfenwick 2200 if ((axis_offset<0) || (axis_offset>getDataPointRank()))
1756     {
1757     throw DataException("Error - Data::trace, axis_offset must be between 0 and rank-2 inclusive.");
1758     }
1759 jfenwick 1796 DataTypes::ShapeType s=getDataPointShape();
1760 ksteube 775 if (getDataPointRank()==2) {
1761 jfenwick 1796 DataTypes::ShapeType ev_shape;
1762 ksteube 775 Data ev(0.,ev_shape,getFunctionSpace());
1763     ev.typeMatchRight(*this);
1764 gross 800 m_data->trace(ev.m_data.get(), axis_offset);
1765 ksteube 775 return ev;
1766     }
1767     if (getDataPointRank()==3) {
1768 jfenwick 1796 DataTypes::ShapeType ev_shape;
1769 ksteube 775 if (axis_offset==0) {
1770     int s2=s[2];
1771     ev_shape.push_back(s2);
1772     }
1773     else if (axis_offset==1) {
1774     int s0=s[0];
1775     ev_shape.push_back(s0);
1776     }
1777     Data ev(0.,ev_shape,getFunctionSpace());
1778     ev.typeMatchRight(*this);
1779 gross 800 m_data->trace(ev.m_data.get(), axis_offset);
1780 ksteube 775 return ev;
1781     }
1782     if (getDataPointRank()==4) {
1783 jfenwick 1796 DataTypes::ShapeType ev_shape;
1784 ksteube 775 if (axis_offset==0) {
1785     ev_shape.push_back(s[2]);
1786     ev_shape.push_back(s[3]);
1787     }
1788     else if (axis_offset==1) {
1789     ev_shape.push_back(s[0]);
1790     ev_shape.push_back(s[3]);
1791     }
1792     else if (axis_offset==2) {
1793     ev_shape.push_back(s[0]);
1794     ev_shape.push_back(s[1]);
1795     }
1796     Data ev(0.,ev_shape,getFunctionSpace());
1797     ev.typeMatchRight(*this);
1798 gross 800 m_data->trace(ev.m_data.get(), axis_offset);
1799 ksteube 775 return ev;
1800     }
1801     else {
1802 gross 800 throw DataException("Error - Data::trace can only be calculated for rank 2, 3 or 4 object.");
1803 ksteube 775 }
1804     }
1805    
1806     Data
1807     Data::transpose(int axis_offset) const
1808 jfenwick 2005 {
1809 jfenwick 2146 MAKELAZYOPOFF(TRANS,axis_offset)
1810 jfenwick 1796 DataTypes::ShapeType s=getDataPointShape();
1811     DataTypes::ShapeType ev_shape;
1812 ksteube 775 // Here's the equivalent of python s_out=s[axis_offset:]+s[:axis_offset]
1813     // which goes thru all shape vector elements starting with axis_offset (at index=rank wrap around to 0)
1814     int rank=getDataPointRank();
1815     if (axis_offset<0 || axis_offset>rank) {
1816     throw DataException("Error - Data::transpose must have 0 <= axis_offset <= rank=" + rank);
1817     }
1818     for (int i=0; i<rank; i++) {
1819     int index = (axis_offset+i)%rank;
1820     ev_shape.push_back(s[index]); // Append to new shape
1821     }
1822     Data ev(0.,ev_shape,getFunctionSpace());
1823     ev.typeMatchRight(*this);
1824     m_data->transpose(ev.m_data.get(), axis_offset);
1825     return ev;
1826 jgs 123 }
1827    
1828 gross 576 Data
1829     Data::eigenvalues() const
1830     {
1831 jfenwick 2005 if (isLazy())
1832     {
1833     Data temp(*this); // to get around the fact that you can't resolve a const Data
1834     temp.resolve();
1835     return temp.eigenvalues();
1836     }
1837 gross 576 // check input
1838 jfenwick 1796 DataTypes::ShapeType s=getDataPointShape();
1839 ksteube 1312 if (getDataPointRank()!=2)
1840 gross 576 throw DataException("Error - Data::eigenvalues can only be calculated for rank 2 object.");
1841 ksteube 1312 if(s[0] != s[1])
1842 gross 576 throw DataException("Error - Data::eigenvalues can only be calculated for object with equal first and second dimension.");
1843     // create return
1844 jfenwick 1796 DataTypes::ShapeType ev_shape(1,s[0]);
1845 gross 576 Data ev(0.,ev_shape,getFunctionSpace());
1846     ev.typeMatchRight(*this);
1847     m_data->eigenvalues(ev.m_data.get());
1848     return ev;
1849     }
1850    
1851 jgs 121 const boost::python::tuple
1852 gross 576 Data::eigenvalues_and_eigenvectors(const double tol) const
1853     {
1854 jfenwick 2005 if (isLazy())
1855     {
1856     Data temp(*this); // to get around the fact that you can't resolve a const Data
1857     temp.resolve();
1858     return temp.eigenvalues_and_eigenvectors(tol);
1859     }
1860 jfenwick 1796 DataTypes::ShapeType s=getDataPointShape();
1861 ksteube 1312 if (getDataPointRank()!=2)
1862 gross 576 throw DataException("Error - Data::eigenvalues and eigenvectors can only be calculated for rank 2 object.");
1863 ksteube 1312 if(s[0] != s[1])
1864 gross 576 throw DataException("Error - Data::eigenvalues and eigenvectors can only be calculated for object with equal first and second dimension.");
1865     // create return
1866 jfenwick 1796 DataTypes::ShapeType ev_shape(1,s[0]);
1867 gross 576 Data ev(0.,ev_shape,getFunctionSpace());
1868     ev.typeMatchRight(*this);
1869 jfenwick 1796 DataTypes::ShapeType V_shape(2,s[0]);
1870 gross 576 Data V(0.,V_shape,getFunctionSpace());
1871     V.typeMatchRight(*this);
1872     m_data->eigenvalues_and_eigenvectors(ev.m_data.get(),V.m_data.get(),tol);
1873     return make_tuple(boost::python::object(ev),boost::python::object(V));
1874     }
1875    
1876     const boost::python::tuple
1877 gross 921 Data::minGlobalDataPoint() const
1878 jgs 121 {
1879 gross 921 // NB: calc_minGlobalDataPoint( had to be split off from minGlobalDataPoint( as boost::make_tuple causes an
1880 jgs 148 // abort (for unknown reasons) if there are openmp directives with it in the
1881     // surrounding function
1882    
1883     int DataPointNo;
1884 gross 921 int ProcNo;
1885     calc_minGlobalDataPoint(ProcNo,DataPointNo);
1886     return make_tuple(ProcNo,DataPointNo);
1887 jgs 148 }
1888    
1889     void
1890 gross 921 Data::calc_minGlobalDataPoint(int& ProcNo,
1891     int& DataPointNo) const
1892 jgs 148 {
1893 jfenwick 2005 if (isLazy())
1894     {
1895     Data temp(*this); // to get around the fact that you can't resolve a const Data
1896     temp.resolve();
1897     return temp.calc_minGlobalDataPoint(ProcNo,DataPointNo);
1898     }
1899 jgs 148 int i,j;
1900     int lowi=0,lowj=0;
1901     double min=numeric_limits<double>::max();
1902    
1903 jgs 121 Data temp=minval();
1904    
1905     int numSamples=temp.getNumSamples();
1906     int numDPPSample=temp.getNumDataPointsPerSample();
1907    
1908 jgs 148 double next,local_min;
1909 jfenwick 1946 int local_lowi=0,local_lowj=0;
1910 jgs 121
1911 caltinay 2081 #pragma omp parallel firstprivate(local_lowi,local_lowj) private(next,local_min)
1912 jgs 148 {
1913     local_min=min;
1914     #pragma omp for private(i,j) schedule(static)
1915     for (i=0; i<numSamples; i++) {
1916     for (j=0; j<numDPPSample; j++) {
1917 jfenwick 1796 next=temp.getDataAtOffset(temp.getDataOffset(i,j));
1918 jgs 148 if (next<local_min) {
1919     local_min=next;
1920     local_lowi=i;
1921     local_lowj=j;
1922     }
1923 jgs 121 }
1924     }
1925 jgs 148 #pragma omp critical
1926     if (local_min<min) {
1927     min=local_min;
1928     lowi=local_lowi;
1929     lowj=local_lowj;
1930     }
1931 jgs 121 }
1932    
1933 bcumming 782 #ifdef PASO_MPI
1934     // determine the processor on which the minimum occurs
1935 jfenwick 1796 next = temp.getDataPoint(lowi,lowj);
1936 bcumming 782 int lowProc = 0;
1937     double *globalMins = new double[get_MPISize()+1];
1938 caltinay 2081 int error;
1939     error = MPI_Gather ( &next, 1, MPI_DOUBLE, globalMins, 1, MPI_DOUBLE, 0, get_MPIComm() );
1940 ksteube 1312
1941 bcumming 782 if( get_MPIRank()==0 ){
1942     next = globalMins[lowProc];
1943     for( i=1; i<get_MPISize(); i++ )
1944     if( next>globalMins[i] ){
1945     lowProc = i;
1946     next = globalMins[i];
1947     }
1948     }
1949     MPI_Bcast( &lowProc, 1, MPI_DOUBLE, 0, get_MPIComm() );
1950    
1951     delete [] globalMins;
1952     ProcNo = lowProc;
1953 bcumming 790 #else
1954     ProcNo = 0;
1955 bcumming 782 #endif
1956 gross 921 DataPointNo = lowj + lowi * numDPPSample;
1957 jgs 121 }
1958    
1959 jgs 104 void
1960     Data::saveDX(std::string fileName) const
1961     {
1962 jfenwick 1803 if (isEmpty())
1963     {
1964     throw DataException("Error - Operations not permitted on instances of DataEmpty.");
1965     }
1966 jfenwick 2005 if (isLazy())
1967     {
1968     Data temp(*this); // to get around the fact that you can't resolve a const Data
1969     temp.resolve();
1970     temp.saveDX(fileName);
1971     return;
1972     }
1973 jgs 153 boost::python::dict args;
1974     args["data"]=boost::python::object(this);
1975 jfenwick 1872 getDomain()->saveDX(fileName,args);
1976 jgs 104 return;
1977     }
1978    
1979 jgs 110 void
1980     Data::saveVTK(std::string fileName) const
1981     {
1982 jfenwick 1803 if (isEmpty())
1983     {
1984     throw DataException("Error - Operations not permitted on instances of DataEmpty.");
1985     }
1986 jfenwick 2005 if (isLazy())
1987     {
1988     Data temp(*this); // to get around the fact that you can't resolve a const Data
1989     temp.resolve();
1990     temp.saveVTK(fileName);
1991     return;
1992     }
1993 jgs 153 boost::python::dict args;
1994     args["data"]=boost::python::object(this);
1995 jfenwick 1872 getDomain()->saveVTK(fileName,args);
1996 jgs 110 return;
1997     }
1998    
1999 jgs 102 Data&
2000     Data::operator+=(const Data& right)
2001     {
2002 gross 783 if (isProtected()) {
2003     throw DataException("Error - attempt to update protected Data object.");
2004     }
2005 jfenwick 2146 MAKELAZYBINSELF(right,ADD) // for lazy + is equivalent to +=
2006     binaryOp(right,plus<double>());
2007     return (*this);
2008 jgs 94 }
2009    
2010 jgs 102 Data&
2011     Data::operator+=(const boost::python::object& right)
2012 jgs 94 {
2013 jfenwick 2005 if (isProtected()) {
2014     throw DataException("Error - attempt to update protected Data object.");
2015     }
2016 gross 854 Data tmp(right,getFunctionSpace(),false);
2017 jfenwick 2146 MAKELAZYBINSELF(tmp,ADD)
2018     binaryOp(tmp,plus<double>());
2019     return (*this);
2020 jgs 94 }
2021 jfenwick 2005
2022     // Hmmm, operator= makes a deep copy but the copy constructor does not?
2023 ksteube 1312 Data&
2024     Data::operator=(const Data& other)
2025     {
2026     copy(other);
2027     return (*this);
2028     }
2029 jgs 94
2030 jgs 102 Data&
2031     Data::operator-=(const Data& right)
2032 jgs 94 {
2033 gross 783 if (isProtected()) {
2034     throw DataException("Error - attempt to update protected Data object.");
2035     }
2036 jfenwick 2146 MAKELAZYBINSELF(right,SUB)
2037     binaryOp(right,minus<double>());
2038     return (*this);
2039 jgs 94 }
2040    
2041 jgs 102 Data&
2042     Data::operator-=(const boost::python::object& right)
2043 jgs 94 {
2044 jfenwick 2005 if (isProtected()) {
2045     throw DataException("Error - attempt to update protected Data object.");
2046     }
2047 gross 854 Data tmp(right,getFunctionSpace(),false);
2048 jfenwick 2146 MAKELAZYBINSELF(tmp,SUB)
2049     binaryOp(tmp,minus<double>());
2050     return (*this);
2051 jgs 94 }
2052    
2053 jgs 102 Data&
2054     Data::operator*=(const Data& right)
2055 jgs 94 {
2056 gross 783 if (isProtected()) {
2057     throw DataException("Error - attempt to update protected Data object.");
2058     }
2059 jfenwick 2146 MAKELAZYBINSELF(right,MUL)
2060     binaryOp(right,multiplies<double>());
2061     return (*this);
2062 jgs 94 }
2063    
2064 jgs 102 Data&
2065     Data::operator*=(const boost::python::object& right)
2066 jfenwick 2005 {
2067     if (isProtected()) {
2068     throw DataException("Error - attempt to update protected Data object.");
2069     }
2070 gross 854 Data tmp(right,getFunctionSpace(),false);
2071 jfenwick 2146 MAKELAZYBINSELF(tmp,MUL)
2072     binaryOp(tmp,multiplies<double>());
2073     return (*this);
2074 jgs 94 }
2075    
2076 jgs 102 Data&
2077     Data::operator/=(const Data& right)
2078 jgs 94 {
2079 gross 783 if (isProtected()) {
2080     throw DataException("Error - attempt to update protected Data object.");
2081     }
2082 jfenwick 2146 MAKELAZYBINSELF(right,DIV)
2083     binaryOp(right,divides<double>());
2084     return (*this);
2085 jgs 94 }
2086    
2087 jgs 102 Data&
2088     Data::operator/=(const boost::python::object& right)
2089 jgs 94 {
2090 jfenwick 2005 if (isProtected()) {
2091     throw DataException("Error - attempt to update protected Data object.");
2092     }
2093 gross 854 Data tmp(right,getFunctionSpace(),false);
2094 jfenwick 2146 MAKELAZYBINSELF(tmp,DIV)
2095     binaryOp(tmp,divides<double>());
2096     return (*this);
2097 jgs 94 }
2098    
2099 jgs 102 Data
2100 gross 699 Data::rpowO(const boost::python::object& left) const
2101     {
2102     Data left_d(left,*this);
2103     return left_d.powD(*this);
2104     }
2105    
2106     Data
2107 jgs 102 Data::powO(const boost::python::object& right) const
2108 jgs 94 {
2109 gross 854 Data tmp(right,getFunctionSpace(),false);
2110     return powD(tmp);
2111 jgs 94 }
2112    
2113 jgs 102 Data
2114     Data::powD(const Data& right) const
2115 jgs 94 {
2116 jfenwick 2146 MAKELAZYBIN(right,POW)
2117 matt 1350 return C_TensorBinaryOperation<double (*)(double, double)>(*this, right, ::pow);
2118 jgs 94 }
2119    
2120     //
2121 jgs 123 // NOTE: It is essential to specify the namespace this operator belongs to
2122 jgs 102 Data
2123     escript::operator+(const Data& left, const Data& right)
2124 jgs 94 {
2125 jfenwick 2146 MAKELAZYBIN2(left,right,ADD)
2126 matt 1327 return C_TensorBinaryOperation(left, right, plus<double>());
2127 jgs 94 }
2128    
2129     //
2130 jgs 123 // NOTE: It is essential to specify the namespace this operator belongs to
2131 jgs 102 Data
2132     escript::operator-(const Data& left, const Data& right)
2133 jgs 94 {
2134 jfenwick 2146 MAKELAZYBIN2(left,right,SUB)
2135 matt 1327 return C_TensorBinaryOperation(left, right, minus<double>());
2136 jgs 94 }
2137    
2138     //
2139 jgs 123 // NOTE: It is essential to specify the namespace this operator belongs to
2140 jgs 102 Data
2141     escript::operator*(const Data& left, const Data& right)
2142 jgs 94 {
2143 jfenwick 2146 MAKELAZYBIN2(left,right,MUL)
2144 matt 1327 return C_TensorBinaryOperation(left, right, multiplies<double>());
2145 jgs 94 }
2146    
2147     //
2148 jgs 123 // NOTE: It is essential to specify the namespace this operator belongs to
2149 jgs 102 Data
2150     escript::operator/(const Data& left, const Data& right)
2151 jgs 94 {
2152 jfenwick 2146 MAKELAZYBIN2(left,right,DIV)
2153 matt 1327 return C_TensorBinaryOperation(left, right, divides<double>());
2154 jgs 94 }
2155    
2156     //
2157 jgs 123 // NOTE: It is essential to specify the namespace this operator belongs to
2158 jgs 102 Data
2159     escript::operator+(const Data& left, const boost::python::object& right)
2160 jgs 94 {
2161 jfenwick 2146 Data tmp(right,left.getFunctionSpace(),false);
2162     MAKELAZYBIN2(left,tmp,ADD)
2163     return left+tmp;
2164 jgs 94 }
2165    
2166     //
2167 jgs 123 // NOTE: It is essential to specify the namespace this operator belongs to
2168 jgs 102 Data
2169     escript::operator-(const Data& left, const boost::python::object& right)
2170 jgs 94 {
2171 jfenwick 2146 Data tmp(right,left.getFunctionSpace(),false);
2172     MAKELAZYBIN2(left,tmp,SUB)
2173     return left-tmp;
2174 jgs 94 }
2175    
2176     //
2177 jgs 123 // NOTE: It is essential to specify the namespace this operator belongs to
2178 jgs 102 Data
2179     escript::operator*(const Data& left, const boost::python::object& right)
2180 jgs 94 {
2181 jfenwick 2146 Data tmp(right,left.getFunctionSpace(),false);
2182     MAKELAZYBIN2(left,tmp,MUL)
2183     return left*tmp;
2184 jgs 94 }
2185    
2186     //
2187 jgs 123 // NOTE: It is essential to specify the namespace this operator belongs to
2188 jgs 102 Data
2189     escript::operator/(const Data& left, const boost::python::object& right)
2190 jgs 94 {
2191 jfenwick 2146 Data tmp(right,left.getFunctionSpace(),false);
2192     MAKELAZYBIN2(left,tmp,DIV)
2193     return left/tmp;
2194 jgs 94 }
2195    
2196     //
2197 jgs 123 // NOTE: It is essential to specify the namespace this operator belongs to
2198 jgs 102 Data
2199     escript::operator+(const boost::python::object& left, const Data& right)
2200 jgs 94 {
2201 jfenwick 2146 Data tmp(left,right.getFunctionSpace(),false);
2202     MAKELAZYBIN2(tmp,right,ADD)
2203     return tmp+right;
2204 jgs 94 }
2205    
2206     //
2207 jgs 123 // NOTE: It is essential to specify the namespace this operator belongs to
2208 jgs 102 Data
2209     escript::operator-(const boost::python::object& left, const Data& right)
2210 jgs 94 {
2211 jfenwick 2146 Data tmp(left,right.getFunctionSpace(),false);
2212     MAKELAZYBIN2(tmp,right,SUB)
2213     return tmp-right;
2214 jgs 94 }
2215    
2216     //
2217 jgs 123 // NOTE: It is essential to specify the namespace this operator belongs to
2218 jgs 102 Data
2219     escript::operator*(const boost::python::object& left, const Data& right)
2220 jgs 94 {
2221 jfenwick 2146 Data tmp(left,right.getFunctionSpace(),false);
2222     MAKELAZYBIN2(tmp,right,MUL)
2223     return tmp*right;
2224 jgs 94 }
2225    
2226     //
2227 jgs 123 // NOTE: It is essential to specify the namespace this operator belongs to
2228 jgs 102 Data
2229     escript::operator/(const boost::python::object& left, const Data& right)
2230 jgs 94 {
2231 jfenwick 2146 Data tmp(left,right.getFunctionSpace(),false);
2232     MAKELAZYBIN2(tmp,right,DIV)
2233     return tmp/right;
2234 jgs 94 }
2235    
2236 jgs 102
2237 bcumming 751 /* TODO */
2238     /* global reduction */
2239 jgs 102 Data
2240 ksteube 1312 Data::getItem(const boost::python::object& key) const
2241 jgs 94 {
2242    
2243 jfenwick 1796 DataTypes::RegionType slice_region=DataTypes::getSliceRegion(getDataPointShape(),key);
2244 jgs 94
2245 jfenwick 1796 if (slice_region.size()!=getDataPointRank()) {
2246 jgs 102 throw DataException("Error - slice size does not match Data rank.");
2247 jgs 94 }
2248    
2249 jgs 102 return getSlice(slice_region);
2250 jgs 94 }
2251    
2252 bcumming 751 /* TODO */
2253     /* global reduction */
2254 jgs 94 Data
2255 jfenwick 1796 Data::getSlice(const DataTypes::RegionType& region) const
2256 jgs 94 {
2257 jgs 102 return Data(*this,region);
2258 jgs 94 }
2259    
2260 bcumming 751 /* TODO */
2261     /* global reduction */
2262 jgs 94 void
2263 jgs 102 Data::setItemO(const boost::python::object& key,
2264     const boost::python::object& value)
2265 jgs 94 {
2266 jgs 102 Data tempData(value,getFunctionSpace());
2267     setItemD(key,tempData);
2268     }
2269    
2270     void
2271     Data::setItemD(const boost::python::object& key,
2272     const Data& value)
2273     {
2274 jfenwick 1796 // const DataArrayView& view=getPointDataView();
2275 jgs 123
2276 jfenwick 1796 DataTypes::RegionType slice_region=DataTypes::getSliceRegion(getDataPointShape(),key);
2277     if (slice_region.size()!=getDataPointRank()) {
2278 jgs 94 throw DataException("Error - slice size does not match Data rank.");
2279     }
2280 jgs 108 if (getFunctionSpace()!=value.getFunctionSpace()) {
2281     setSlice(Data(value,getFunctionSpace()),slice_region);
2282     } else {
2283     setSlice(value,slice_region);
2284     }
2285 jgs 94 }
2286    
2287     void
2288 jgs 102 Data::setSlice(const Data& value,
2289 jfenwick 1796 const DataTypes::RegionType& region)
2290 jgs 94 {
2291 gross 783 if (isProtected()) {
2292     throw DataException("Error - attempt to update protected Data object.");
2293     }
2294 jfenwick 2005 FORCERESOLVE;
2295     /* if (isLazy())
2296     {
2297     throw DataException("Error - setSlice not permitted on lazy data.");
2298     }*/
2299 jgs 102 Data tempValue(value);
2300     typeMatchLeft(tempValue);
2301     typeMatchRight(tempValue);
2302 jfenwick 2005 getReady()->setSlice(tempValue.m_data.get(),region);
2303 jgs 102 }
2304    
2305     void
2306     Data::typeMatchLeft(Data& right) const
2307     {
2308 jfenwick 2005 if (right.isLazy() && !isLazy())
2309     {
2310     right.resolve();
2311     }
2312 jgs 102 if (isExpanded()){
2313     right.expand();
2314     } else if (isTagged()) {
2315     if (right.isConstant()) {
2316     right.tag();
2317     }
2318     }
2319     }
2320    
2321     void
2322     Data::typeMatchRight(const Data& right)
2323     {
2324 jfenwick 2005 if (isLazy() && !right.isLazy())
2325     {
2326     resolve();
2327     }
2328 jgs 94 if (isTagged()) {
2329     if (right.isExpanded()) {
2330     expand();
2331     }
2332     } else if (isConstant()) {
2333     if (right.isExpanded()) {
2334     expand();
2335     } else if (right.isTagged()) {
2336     tag();
2337     }
2338     }
2339     }
2340    
2341     void
2342 gross 1044 Data::setTaggedValueByName(std::string name,
2343 ksteube 1312 const boost::python::object& value)
2344 gross 1044 {
2345 jfenwick 1872 if (getFunctionSpace().getDomain()->isValidTagName(name)) {
2346 jfenwick 2005 FORCERESOLVE;
2347 jfenwick 1872 int tagKey=getFunctionSpace().getDomain()->getTag(name);
2348 gross 1044 setTaggedValue(tagKey,value);
2349     }
2350     }
2351     void
2352 jgs 94 Data::setTaggedValue(int tagKey,
2353     const boost::python::object& value)
2354     {
2355 gross 783 if (isProtected()) {
2356     throw DataException("Error - attempt to update protected Data object.");
2357     }
2358 jgs 94 //
2359     // Ensure underlying data object is of type DataTagged
2360 jfenwick 2005 FORCERESOLVE;
2361 gross 1358 if (isConstant()) tag();
2362 matt 1319 numeric::array asNumArray(value);
2363 jgs 94
2364 matt 1319 // extract the shape of the numarray
2365 jfenwick 1796 DataTypes::ShapeType tempShape;
2366 matt 1319 for (int i=0; i < asNumArray.getrank(); i++) {
2367     tempShape.push_back(extract<int>(asNumArray.getshape()[i]));
2368     }
2369    
2370 jfenwick 1796 DataVector temp_data2;
2371 jfenwick 2105 temp_data2.copyFromNumArray(asNumArray,1);
2372 jfenwick 1796
2373 jfenwick 2005 m_data->setTaggedValue(tagKey,tempShape, temp_data2);
2374 jgs 94 }
2375    
2376 jfenwick 1796
2377 jgs 110 void
2378 jgs 121 Data::setTaggedValueFromCPP(int tagKey,
2379 jfenwick 1796 const DataTypes::ShapeType& pointshape,
2380     const DataTypes::ValueType& value,
2381     int dataOffset)
2382 jgs 121 {
2383 gross 783 if (isProtected()) {
2384     throw DataException("Error - attempt to update protected Data object.");
2385     }
2386 jgs 121 //
2387     // Ensure underlying data object is of type DataTagged
2388 jfenwick 2005 FORCERESOLVE;
2389 gross 1358 if (isConstant()) tag();
2390 jgs 121 //
2391     // Call DataAbstract::setTaggedValue
2392 jfenwick 1796 m_data->setTaggedValue(tagKey,pointshape, value, dataOffset);
2393 jgs 121 }
2394    
2395 jgs 149 int
2396     Data::getTagNumber(int dpno)
2397     {
2398 jfenwick 1803 if (isEmpty())
2399     {
2400     throw DataException("Error - operation not permitted on instances of DataEmpty.");
2401     }
2402 gross 1409 return getFunctionSpace().getTagFromDataPointNo(dpno);
2403 jgs 149 }
2404    
2405 jgs 119
2406 jgs 94 ostream& escript::operator<<(ostream& o, const Data& data)
2407     {
2408     o << data.toString();
2409     return o;
2410     }
2411 bcumming 782
2412 ksteube 813 Data
2413     escript::C_GeneralTensorProduct(Data& arg_0,
2414     Data& arg_1,
2415     int axis_offset,
2416     int transpose)
2417     {
2418 gross 826 // General tensor product: res(SL x SR) = arg_0(SL x SM) * arg_1(SM x SR)
2419 ksteube 813 // SM is the product of the last axis_offset entries in arg_0.getShape().
2420    
2421 jfenwick 2005 // deal with any lazy data
2422 jfenwick 2066 // if (arg_0.isLazy()) {arg_0.resolve();}
2423     // if (arg_1.isLazy()) {arg_1.resolve();}
2424 jfenwick 2199 if (arg_0.isLazy() || arg_1.isLazy() || (AUTOLAZYON && (arg_0.isExpanded() || arg_1.isExpanded())))
2425 jfenwick 2066 {
2426     DataLazy* c=new DataLazy(arg_0.borrowDataPtr(), arg_1.borrowDataPtr(), PROD, axis_offset,transpose);
2427     return Data(c);
2428     }
2429 jfenwick 2005
2430 ksteube 813 // Interpolate if necessary and find an appropriate function space
2431 gross 826 Data arg_0_Z, arg_1_Z;
2432 ksteube 813 if (arg_0.getFunctionSpace()!=arg_1.getFunctionSpace()) {
2433     if (arg_0.probeInterpolation(arg_1.getFunctionSpace())) {
2434 gross 826 arg_0_Z = arg_0.interpolate(arg_1.getFunctionSpace());
2435     arg_1_Z = Data(arg_1);
2436 ksteube 813 }
2437     else if (arg_1.probeInterpolation(arg_0.getFunctionSpace())) {
2438 gross 826 arg_1_Z=arg_1.interpolate(arg_0.getFunctionSpace());
2439     arg_0_Z =Data(arg_0);
2440 ksteube 813 }
2441     else {
2442     throw DataException("Error - C_GeneralTensorProduct: arguments have incompatible function spaces.");
2443     }
2444 gross 826 } else {
2445     arg_0_Z = Data(arg_0);
2446     arg_1_Z = Data(arg_1);
2447 ksteube 813 }
2448     // Get rank and shape of inputs
2449 gross 826 int