1 |
|
2 |
/******************************************************* |
3 |
* |
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 |
|
14 |
|
15 |
#include "Data.h" |
16 |
|
17 |
#include "DataExpanded.h" |
18 |
#include "DataConstant.h" |
19 |
#include "DataTagged.h" |
20 |
#include "DataEmpty.h" |
21 |
#include "DataLazy.h" |
22 |
#include "FunctionSpaceFactory.h" |
23 |
#include "AbstractContinuousDomain.h" |
24 |
#include "UnaryFuncs.h" |
25 |
#include "FunctionSpaceException.h" |
26 |
#include "EscriptParams.h" |
27 |
|
28 |
extern "C" { |
29 |
#include "esysUtils/blocktimer.h" |
30 |
} |
31 |
|
32 |
#include <fstream> |
33 |
#include <algorithm> |
34 |
#include <vector> |
35 |
#include <functional> |
36 |
#include <sstream> // so we can throw messages about ranks |
37 |
|
38 |
#include <boost/python/dict.hpp> |
39 |
#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 |
// 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 |
|
51 |
Data::Data() |
52 |
{ |
53 |
// |
54 |
// Default data is type DataEmpty |
55 |
DataAbstract* temp=new DataEmpty(); |
56 |
m_data=temp->getPtr(); |
57 |
m_protected=false; |
58 |
} |
59 |
|
60 |
Data::Data(double value, |
61 |
const tuple& shape, |
62 |
const FunctionSpace& what, |
63 |
bool expanded) |
64 |
{ |
65 |
DataTypes::ShapeType dataPointShape; |
66 |
for (int i = 0; i < shape.attr("__len__")(); ++i) { |
67 |
dataPointShape.push_back(extract<const int>(shape[i])); |
68 |
} |
69 |
|
70 |
int len = DataTypes::noValues(dataPointShape); |
71 |
DataVector temp_data(len,value,len); |
72 |
initialise(temp_data, dataPointShape, what, expanded); |
73 |
m_protected=false; |
74 |
} |
75 |
|
76 |
Data::Data(double value, |
77 |
const DataTypes::ShapeType& dataPointShape, |
78 |
const FunctionSpace& what, |
79 |
bool expanded) |
80 |
{ |
81 |
int len = DataTypes::noValues(dataPointShape); |
82 |
|
83 |
DataVector temp_data(len,value,len); |
84 |
// DataArrayView temp_dataView(temp_data, dataPointShape); |
85 |
|
86 |
// initialise(temp_dataView, what, expanded); |
87 |
initialise(temp_data, dataPointShape, what, expanded); |
88 |
|
89 |
m_protected=false; |
90 |
} |
91 |
|
92 |
Data::Data(const Data& inData) |
93 |
{ |
94 |
m_data=inData.m_data; |
95 |
m_protected=inData.isProtected(); |
96 |
} |
97 |
|
98 |
|
99 |
Data::Data(const Data& inData, |
100 |
const DataTypes::RegionType& region) |
101 |
{ |
102 |
DataAbstract_ptr dat=inData.m_data; |
103 |
if (inData.isLazy()) |
104 |
{ |
105 |
dat=inData.m_data->resolve(); |
106 |
} |
107 |
else |
108 |
{ |
109 |
dat=inData.m_data; |
110 |
} |
111 |
// |
112 |
// Create Data which is a slice of another Data |
113 |
DataAbstract* tmp = dat->getSlice(region); |
114 |
m_data=DataAbstract_ptr(tmp); |
115 |
m_protected=false; |
116 |
} |
117 |
|
118 |
Data::Data(const Data& inData, |
119 |
const FunctionSpace& functionspace) |
120 |
{ |
121 |
if (inData.isEmpty()) |
122 |
{ |
123 |
throw DataException("Error - will not interpolate for instances of DataEmpty."); |
124 |
} |
125 |
if (inData.getFunctionSpace()==functionspace) { |
126 |
m_data=inData.m_data; |
127 |
} |
128 |
else |
129 |
{ |
130 |
|
131 |
if (inData.isConstant()) { // for a constant function, we just need to use the new function space |
132 |
if (!inData.probeInterpolation(functionspace)) |
133 |
{ // Even though this is constant, we still need to check whether interpolation is allowed |
134 |
throw FunctionSpaceException("Cannot interpolate across to the domain of the specified FunctionSpace. (DataConstant)"); |
135 |
} |
136 |
// if the data is not lazy, this will just be a cast to DataReady |
137 |
DataReady_ptr dr=inData.m_data->resolve(); |
138 |
DataConstant* dc=new DataConstant(functionspace,inData.m_data->getShape(),dr->getVector()); |
139 |
m_data=DataAbstract_ptr(dc); |
140 |
} else { |
141 |
Data tmp(0,inData.getDataPointShape(),functionspace,true); |
142 |
// Note: Must use a reference or pointer to a derived object |
143 |
// in order to get polymorphic behaviour. Shouldn't really |
144 |
// be able to create an instance of AbstractDomain but that was done |
145 |
// as a boost:python work around which may no longer be required. |
146 |
/*const AbstractDomain& inDataDomain=inData.getDomain();*/ |
147 |
const_Domain_ptr inDataDomain=inData.getDomain(); |
148 |
if (inDataDomain==functionspace.getDomain()) { |
149 |
inDataDomain->interpolateOnDomain(tmp,inData); |
150 |
} else { |
151 |
inDataDomain->interpolateACross(tmp,inData); |
152 |
} |
153 |
m_data=tmp.m_data; |
154 |
} |
155 |
} |
156 |
m_protected=false; |
157 |
} |
158 |
|
159 |
Data::Data(DataAbstract* underlyingdata) |
160 |
{ |
161 |
// m_data=shared_ptr<DataAbstract>(underlyingdata); |
162 |
m_data=underlyingdata->getPtr(); |
163 |
m_protected=false; |
164 |
} |
165 |
|
166 |
Data::Data(DataAbstract_ptr underlyingdata) |
167 |
{ |
168 |
m_data=underlyingdata; |
169 |
m_protected=false; |
170 |
} |
171 |
|
172 |
|
173 |
Data::Data(const numeric::array& value, |
174 |
const FunctionSpace& what, |
175 |
bool expanded) |
176 |
{ |
177 |
initialise(value,what,expanded); |
178 |
m_protected=false; |
179 |
} |
180 |
|
181 |
Data::Data(const DataTypes::ValueType& value, |
182 |
const DataTypes::ShapeType& shape, |
183 |
const FunctionSpace& what, |
184 |
bool expanded) |
185 |
{ |
186 |
initialise(value,shape,what,expanded); |
187 |
m_protected=false; |
188 |
} |
189 |
|
190 |
|
191 |
Data::Data(const object& value, |
192 |
const FunctionSpace& what, |
193 |
bool expanded) |
194 |
{ |
195 |
numeric::array asNumArray(value); |
196 |
initialise(asNumArray,what,expanded); |
197 |
m_protected=false; |
198 |
} |
199 |
|
200 |
|
201 |
Data::Data(const object& value, |
202 |
const Data& other) |
203 |
{ |
204 |
numeric::array asNumArray(value); |
205 |
|
206 |
// extract the shape of the numarray |
207 |
DataTypes::ShapeType tempShape=DataTypes::shapeFromNumArray(asNumArray); |
208 |
if (DataTypes::getRank(tempShape)==0) { |
209 |
|
210 |
|
211 |
// get the space for the data vector |
212 |
int len1 = DataTypes::noValues(tempShape); |
213 |
DataVector temp_data(len1, 0.0, len1); |
214 |
temp_data.copyFromNumArray(asNumArray); |
215 |
|
216 |
int len = DataTypes::noValues(other.getDataPointShape()); |
217 |
|
218 |
DataVector temp2_data(len, temp_data[0]/*temp_dataView()*/, len); |
219 |
DataConstant* t=new DataConstant(other.getFunctionSpace(),other.getDataPointShape(),temp2_data); |
220 |
m_data=DataAbstract_ptr(t); |
221 |
|
222 |
} else { |
223 |
// |
224 |
// Create a DataConstant with the same sample shape as other |
225 |
DataConstant* t=new DataConstant(asNumArray,other.getFunctionSpace()); |
226 |
// boost::shared_ptr<DataAbstract> sp(t); |
227 |
// m_data=sp; |
228 |
m_data=DataAbstract_ptr(t); |
229 |
} |
230 |
m_protected=false; |
231 |
} |
232 |
|
233 |
Data::~Data() |
234 |
{ |
235 |
|
236 |
} |
237 |
|
238 |
|
239 |
|
240 |
void |
241 |
Data::initialise(const boost::python::numeric::array& value, |
242 |
const FunctionSpace& what, |
243 |
bool expanded) |
244 |
{ |
245 |
// |
246 |
// Construct a Data object of the appropriate type. |
247 |
// Construct the object first as there seems to be a bug which causes |
248 |
// undefined behaviour if an exception is thrown during construction |
249 |
// within the shared_ptr constructor. |
250 |
if (expanded) { |
251 |
DataAbstract* temp=new DataExpanded(value, what); |
252 |
// boost::shared_ptr<DataAbstract> temp_data(temp); |
253 |
// m_data=temp_data; |
254 |
m_data=temp->getPtr(); |
255 |
} else { |
256 |
DataAbstract* temp=new DataConstant(value, what); |
257 |
// boost::shared_ptr<DataAbstract> temp_data(temp); |
258 |
// m_data=temp_data; |
259 |
m_data=temp->getPtr(); |
260 |
} |
261 |
} |
262 |
|
263 |
|
264 |
void |
265 |
Data::initialise(const DataTypes::ValueType& value, |
266 |
const DataTypes::ShapeType& shape, |
267 |
const FunctionSpace& what, |
268 |
bool expanded) |
269 |
{ |
270 |
// |
271 |
// Construct a Data object of the appropriate type. |
272 |
// Construct the object first as there seems to be a bug which causes |
273 |
// undefined behaviour if an exception is thrown during construction |
274 |
// within the shared_ptr constructor. |
275 |
if (expanded) { |
276 |
DataAbstract* temp=new DataExpanded(what, shape, value); |
277 |
// boost::shared_ptr<DataAbstract> temp_data(temp); |
278 |
// m_data=temp_data; |
279 |
m_data=temp->getPtr(); |
280 |
} else { |
281 |
DataAbstract* temp=new DataConstant(what, shape, value); |
282 |
// boost::shared_ptr<DataAbstract> temp_data(temp); |
283 |
// m_data=temp_data; |
284 |
m_data=temp->getPtr(); |
285 |
} |
286 |
} |
287 |
|
288 |
|
289 |
// void |
290 |
// Data::CompareDebug(const Data& rd) |
291 |
// { |
292 |
// using namespace std; |
293 |
// bool mismatch=false; |
294 |
// std::cout << "Comparing left and right" << endl; |
295 |
// const DataTagged* left=dynamic_cast<DataTagged*>(m_data.get()); |
296 |
// const DataTagged* right=dynamic_cast<DataTagged*>(rd.m_data.get()); |
297 |
// |
298 |
// if (left==0) |
299 |
// { |
300 |
// cout << "left arg is not a DataTagged\n"; |
301 |
// return; |
302 |
// } |
303 |
// |
304 |
// if (right==0) |
305 |
// { |
306 |
// cout << "right arg is not a DataTagged\n"; |
307 |
// return; |
308 |
// } |
309 |
// cout << "Num elements=" << left->getVector().size() << ":" << right->getVector().size() << std::endl; |
310 |
// cout << "Shapes "; |
311 |
// if (left->getShape()==right->getShape()) |
312 |
// { |
313 |
// cout << "ok\n"; |
314 |
// } |
315 |
// else |
316 |
// { |
317 |
// cout << "Problem: shapes do not match\n"; |
318 |
// mismatch=true; |
319 |
// } |
320 |
// int lim=left->getVector().size(); |
321 |
// if (right->getVector().size()) lim=right->getVector().size(); |
322 |
// for (int i=0;i<lim;++i) |
323 |
// { |
324 |
// if (left->getVector()[i]!=right->getVector()[i]) |
325 |
// { |
326 |
// cout << "[" << i << "] value mismatch " << left->getVector()[i] << ":" << right->getVector()[i] << endl; |
327 |
// mismatch=true; |
328 |
// } |
329 |
// } |
330 |
// |
331 |
// // still need to check the tag map |
332 |
// // also need to watch what is happening to function spaces, are they copied or what? |
333 |
// |
334 |
// const DataTagged::DataMapType& mapleft=left->getTagLookup(); |
335 |
// const DataTagged::DataMapType& mapright=right->getTagLookup(); |
336 |
// |
337 |
// if (mapleft.size()!=mapright.size()) |
338 |
// { |
339 |
// cout << "Maps are different sizes " << mapleft.size() << ":" << mapright.size() << endl; |
340 |
// mismatch=true; |
341 |
// cout << "Left map\n"; |
342 |
// DataTagged::DataMapType::const_iterator i,j; |
343 |
// for (i=mapleft.begin();i!=mapleft.end();++i) { |
344 |
// cout << "(" << i->first << "=>" << i->second << ")\n"; |
345 |
// } |
346 |
// cout << "Right map\n"; |
347 |
// for (i=mapright.begin();i!=mapright.end();++i) { |
348 |
// cout << "(" << i->first << "=>" << i->second << ")\n"; |
349 |
// } |
350 |
// cout << "End map\n"; |
351 |
// |
352 |
// } |
353 |
// |
354 |
// DataTagged::DataMapType::const_iterator i,j; |
355 |
// for (i=mapleft.begin(),j=mapright.begin();i!=mapleft.end() && j!=mapright.end();++i,++j) { |
356 |
// if ((i->first!=j->first) || (i->second!=j->second)) |
357 |
// { |
358 |
// cout << "(" << i->first << "=>" << i->second << ")"; |
359 |
// cout << ":(" << j->first << "=>" << j->second << ") "; |
360 |
// mismatch=true; |
361 |
// } |
362 |
// } |
363 |
// if (mismatch) |
364 |
// { |
365 |
// cout << "#Mismatch\n"; |
366 |
// } |
367 |
// } |
368 |
|
369 |
escriptDataC |
370 |
Data::getDataC() |
371 |
{ |
372 |
escriptDataC temp; |
373 |
temp.m_dataPtr=(void*)this; |
374 |
return temp; |
375 |
} |
376 |
|
377 |
escriptDataC |
378 |
Data::getDataC() const |
379 |
{ |
380 |
escriptDataC temp; |
381 |
temp.m_dataPtr=(void*)this; |
382 |
return temp; |
383 |
} |
384 |
|
385 |
size_t |
386 |
Data::getSampleBufferSize() const |
387 |
{ |
388 |
return m_data->getSampleBufferSize(); |
389 |
} |
390 |
|
391 |
|
392 |
const boost::python::tuple |
393 |
Data::getShapeTuple() const |
394 |
{ |
395 |
const DataTypes::ShapeType& shape=getDataPointShape(); |
396 |
switch(getDataPointRank()) { |
397 |
case 0: |
398 |
return make_tuple(); |
399 |
case 1: |
400 |
return make_tuple(long_(shape[0])); |
401 |
case 2: |
402 |
return make_tuple(long_(shape[0]),long_(shape[1])); |
403 |
case 3: |
404 |
return make_tuple(long_(shape[0]),long_(shape[1]),long_(shape[2])); |
405 |
case 4: |
406 |
return make_tuple(long_(shape[0]),long_(shape[1]),long_(shape[2]),long_(shape[3])); |
407 |
default: |
408 |
throw DataException("Error - illegal Data rank."); |
409 |
} |
410 |
} |
411 |
|
412 |
|
413 |
// The different name is needed because boost has trouble with overloaded functions. |
414 |
// It can't work out what type the function is based soley on its name. |
415 |
// There are ways to fix this involving creating function pointer variables for each form |
416 |
// but there doesn't seem to be a need given that the methods have the same name from the python point of view |
417 |
Data* |
418 |
Data::copySelf() |
419 |
{ |
420 |
DataAbstract* temp=m_data->deepCopy(); |
421 |
return new Data(temp); |
422 |
} |
423 |
|
424 |
void |
425 |
Data::copy(const Data& other) |
426 |
{ |
427 |
DataAbstract* temp=other.m_data->deepCopy(); |
428 |
DataAbstract_ptr p=temp->getPtr(); |
429 |
m_data=p; |
430 |
} |
431 |
|
432 |
|
433 |
Data |
434 |
Data::delay() |
435 |
{ |
436 |
DataLazy* dl=new DataLazy(m_data); |
437 |
return Data(dl); |
438 |
} |
439 |
|
440 |
void |
441 |
Data::delaySelf() |
442 |
{ |
443 |
if (!isLazy()) |
444 |
{ |
445 |
m_data=(new DataLazy(m_data))->getPtr(); |
446 |
} |
447 |
} |
448 |
|
449 |
void |
450 |
Data::setToZero() |
451 |
{ |
452 |
if (isEmpty()) |
453 |
{ |
454 |
throw DataException("Error - Operations not permitted on instances of DataEmpty."); |
455 |
} |
456 |
m_data->setToZero(); |
457 |
} |
458 |
|
459 |
void |
460 |
Data::copyWithMask(const Data& other, |
461 |
const Data& mask) |
462 |
{ |
463 |
// 1. Interpolate if required so all Datas use the same FS as this |
464 |
// 2. Tag or Expand so that all Data's are the same type |
465 |
// 3. Iterate over the data vectors copying values where mask is >0 |
466 |
if (other.isEmpty() || mask.isEmpty()) |
467 |
{ |
468 |
throw DataException("Error - copyWithMask not permitted using instances of DataEmpty."); |
469 |
} |
470 |
Data other2(other); |
471 |
Data mask2(mask); |
472 |
other2.resolve(); |
473 |
mask2.resolve(); |
474 |
this->resolve(); |
475 |
FunctionSpace myFS=getFunctionSpace(); |
476 |
FunctionSpace oFS=other2.getFunctionSpace(); |
477 |
FunctionSpace mFS=mask2.getFunctionSpace(); |
478 |
if (oFS!=myFS) |
479 |
{ |
480 |
if (other2.probeInterpolation(myFS)) |
481 |
{ |
482 |
other2=other2.interpolate(myFS); |
483 |
} |
484 |
else |
485 |
{ |
486 |
throw DataException("Error - copyWithMask: other FunctionSpace is not compatible with this one."); |
487 |
} |
488 |
} |
489 |
if (mFS!=myFS) |
490 |
{ |
491 |
if (mask2.probeInterpolation(myFS)) |
492 |
{ |
493 |
mask2=mask2.interpolate(myFS); |
494 |
} |
495 |
else |
496 |
{ |
497 |
throw DataException("Error - copyWithMask: mask FunctionSpace is not compatible with this one."); |
498 |
} |
499 |
} |
500 |
// Ensure that all args have the same type |
501 |
if (this->isExpanded() || mask2.isExpanded() || other2.isExpanded()) |
502 |
{ |
503 |
this->expand(); |
504 |
other2.expand(); |
505 |
mask2.expand(); |
506 |
} |
507 |
else if (this->isTagged() || mask2.isTagged() || other2.isTagged()) |
508 |
{ |
509 |
this->tag(); |
510 |
other2.tag(); |
511 |
mask2.tag(); |
512 |
} |
513 |
else if (this->isConstant() && mask2.isConstant() && other2.isConstant()) |
514 |
{ |
515 |
} |
516 |
else |
517 |
{ |
518 |
throw DataException("Error - Unknown DataAbstract passed to copyWithMask."); |
519 |
} |
520 |
// Now we iterate over the elements |
521 |
DataVector& self=getReadyPtr()->getVector(); |
522 |
const DataVector& ovec=other2.getReadyPtr()->getVector(); |
523 |
const DataVector& mvec=mask2.getReadyPtr()->getVector(); |
524 |
if ((self.size()!=ovec.size()) || (self.size()!=mvec.size())) |
525 |
{ |
526 |
throw DataException("Error - size mismatch in arguments to copyWithMask."); |
527 |
} |
528 |
size_t num_points=self.size(); |
529 |
|
530 |
// OPENMP 3.0 allows unsigned loop vars. |
531 |
#if defined(_OPENMP) && (_OPENMP < 200805) |
532 |
long i; |
533 |
#else |
534 |
size_t i; |
535 |
#endif |
536 |
#pragma omp parallel for private(i) schedule(static) |
537 |
for (i=0;i<num_points;++i) |
538 |
{ |
539 |
if (mvec[i]>0) |
540 |
{ |
541 |
self[i]=ovec[i]; |
542 |
} |
543 |
} |
544 |
} |
545 |
|
546 |
|
547 |
|
548 |
bool |
549 |
Data::isExpanded() const |
550 |
{ |
551 |
DataExpanded* temp=dynamic_cast<DataExpanded*>(m_data.get()); |
552 |
return (temp!=0); |
553 |
} |
554 |
|
555 |
bool |
556 |
Data::actsExpanded() const |
557 |
{ |
558 |
return m_data->actsExpanded(); |
559 |
|
560 |
} |
561 |
|
562 |
|
563 |
bool |
564 |
Data::isTagged() const |
565 |
{ |
566 |
DataTagged* temp=dynamic_cast<DataTagged*>(m_data.get()); |
567 |
return (temp!=0); |
568 |
} |
569 |
|
570 |
bool |
571 |
Data::isEmpty() const |
572 |
{ |
573 |
DataEmpty* temp=dynamic_cast<DataEmpty*>(m_data.get()); |
574 |
return (temp!=0); |
575 |
} |
576 |
|
577 |
bool |
578 |
Data::isConstant() const |
579 |
{ |
580 |
DataConstant* temp=dynamic_cast<DataConstant*>(m_data.get()); |
581 |
return (temp!=0); |
582 |
} |
583 |
|
584 |
bool |
585 |
Data::isLazy() const |
586 |
{ |
587 |
return m_data->isLazy(); |
588 |
} |
589 |
|
590 |
// at the moment this is synonymous with !isLazy() but that could change |
591 |
bool |
592 |
Data::isReady() const |
593 |
{ |
594 |
return (dynamic_cast<DataReady*>(m_data.get())!=0); |
595 |
} |
596 |
|
597 |
|
598 |
void |
599 |
Data::setProtection() |
600 |
{ |
601 |
m_protected=true; |
602 |
} |
603 |
|
604 |
bool |
605 |
Data::isProtected() const |
606 |
{ |
607 |
return m_protected; |
608 |
} |
609 |
|
610 |
|
611 |
|
612 |
void |
613 |
Data::expand() |
614 |
{ |
615 |
if (isConstant()) { |
616 |
DataConstant* tempDataConst=dynamic_cast<DataConstant*>(m_data.get()); |
617 |
DataAbstract* temp=new DataExpanded(*tempDataConst); |
618 |
// shared_ptr<DataAbstract> temp_data(temp); |
619 |
// m_data=temp_data; |
620 |
m_data=temp->getPtr(); |
621 |
} else if (isTagged()) { |
622 |
DataTagged* tempDataTag=dynamic_cast<DataTagged*>(m_data.get()); |
623 |
DataAbstract* temp=new DataExpanded(*tempDataTag); |
624 |
// shared_ptr<DataAbstract> temp_data(temp); |
625 |
// m_data=temp_data; |
626 |
m_data=temp->getPtr(); |
627 |
} else if (isExpanded()) { |
628 |
// |
629 |
// do nothing |
630 |
} else if (isEmpty()) { |
631 |
throw DataException("Error - Expansion of DataEmpty not possible."); |
632 |
} else if (isLazy()) { |
633 |
resolve(); |
634 |
expand(); // resolve might not give us expanded data |
635 |
} else { |
636 |
throw DataException("Error - Expansion not implemented for this Data type."); |
637 |
} |
638 |
} |
639 |
|
640 |
void |
641 |
Data::tag() |
642 |
{ |
643 |
if (isConstant()) { |
644 |
DataConstant* tempDataConst=dynamic_cast<DataConstant*>(m_data.get()); |
645 |
DataAbstract* temp=new DataTagged(*tempDataConst); |
646 |
// shared_ptr<DataAbstract> temp_data(temp); |
647 |
// m_data=temp_data; |
648 |
m_data=temp->getPtr(); |
649 |
} else if (isTagged()) { |
650 |
// do nothing |
651 |
} else if (isExpanded()) { |
652 |
throw DataException("Error - Creating tag data from DataExpanded not possible."); |
653 |
} else if (isEmpty()) { |
654 |
throw DataException("Error - Creating tag data from DataEmpty not possible."); |
655 |
} else if (isLazy()) { |
656 |
DataAbstract_ptr res=m_data->resolve(); |
657 |
if (m_data->isExpanded()) |
658 |
{ |
659 |
throw DataException("Error - data would resolve to DataExpanded, tagging is not possible."); |
660 |
} |
661 |
m_data=res; |
662 |
tag(); |
663 |
} else { |
664 |
throw DataException("Error - Tagging not implemented for this Data type."); |
665 |
} |
666 |
} |
667 |
|
668 |
void |
669 |
Data::resolve() |
670 |
{ |
671 |
if (isLazy()) |
672 |
{ |
673 |
m_data=m_data->resolve(); |
674 |
} |
675 |
} |
676 |
|
677 |
|
678 |
Data |
679 |
Data::oneOver() const |
680 |
{ |
681 |
if (isLazy()) |
682 |
{ |
683 |
DataLazy* c=new DataLazy(borrowDataPtr(),RECIP); |
684 |
return Data(c); |
685 |
} |
686 |
return C_TensorUnaryOperation(*this, bind1st(divides<double>(),1.)); |
687 |
} |
688 |
|
689 |
Data |
690 |
Data::wherePositive() const |
691 |
{ |
692 |
if (isLazy()) |
693 |
{ |
694 |
DataLazy* c=new DataLazy(borrowDataPtr(),GZ); |
695 |
return Data(c); |
696 |
} |
697 |
return C_TensorUnaryOperation(*this, bind2nd(greater<double>(),0.0)); |
698 |
} |
699 |
|
700 |
Data |
701 |
Data::whereNegative() const |
702 |
{ |
703 |
if (isLazy()) |
704 |
{ |
705 |
DataLazy* c=new DataLazy(borrowDataPtr(),LZ); |
706 |
return Data(c); |
707 |
} |
708 |
return C_TensorUnaryOperation(*this, bind2nd(less<double>(),0.0)); |
709 |
} |
710 |
|
711 |
Data |
712 |
Data::whereNonNegative() const |
713 |
{ |
714 |
if (isLazy()) |
715 |
{ |
716 |
DataLazy* c=new DataLazy(borrowDataPtr(),GEZ); |
717 |
return Data(c); |
718 |
} |
719 |
return C_TensorUnaryOperation(*this, bind2nd(greater_equal<double>(),0.0)); |
720 |
} |
721 |
|
722 |
Data |
723 |
Data::whereNonPositive() const |
724 |
{ |
725 |
if (isLazy()) |
726 |
{ |
727 |
DataLazy* c=new DataLazy(borrowDataPtr(),LEZ); |
728 |
return Data(c); |
729 |
} |
730 |
return C_TensorUnaryOperation(*this, bind2nd(less_equal<double>(),0.0)); |
731 |
} |
732 |
|
733 |
Data |
734 |
Data::whereZero(double tol) const |
735 |
{ |
736 |
Data dataAbs=abs(); |
737 |
return C_TensorUnaryOperation(dataAbs, bind2nd(less_equal<double>(),tol)); |
738 |
} |
739 |
|
740 |
Data |
741 |
Data::whereNonZero(double tol) const |
742 |
{ |
743 |
Data dataAbs=abs(); |
744 |
return C_TensorUnaryOperation(dataAbs, bind2nd(greater<double>(),tol)); |
745 |
} |
746 |
|
747 |
Data |
748 |
Data::interpolate(const FunctionSpace& functionspace) const |
749 |
{ |
750 |
return Data(*this,functionspace); |
751 |
} |
752 |
|
753 |
bool |
754 |
Data::probeInterpolation(const FunctionSpace& functionspace) const |
755 |
{ |
756 |
return getFunctionSpace().probeInterpolation(functionspace); |
757 |
// if (getFunctionSpace()==functionspace) { |
758 |
// return true; |
759 |
// } else { |
760 |
// const_Domain_ptr domain=getDomain(); |
761 |
// if (*domain==*functionspace.getDomain()) { |
762 |
// return domain->probeInterpolationOnDomain(getFunctionSpace().getTypeCode(),functionspace.getTypeCode()); |
763 |
// } else { |
764 |
// return domain->probeInterpolationACross(getFunctionSpace().getTypeCode(),*(functionspace.getDomain()),functionspace.getTypeCode()); |
765 |
// } |
766 |
// } |
767 |
} |
768 |
|
769 |
Data |
770 |
Data::gradOn(const FunctionSpace& functionspace) const |
771 |
{ |
772 |
if (isEmpty()) |
773 |
{ |
774 |
throw DataException("Error - operation not permitted on instances of DataEmpty."); |
775 |
} |
776 |
double blocktimer_start = blocktimer_time(); |
777 |
if (functionspace.getDomain()!=getDomain()) |
778 |
throw DataException("Error - gradient cannot be calculated on different domains."); |
779 |
DataTypes::ShapeType grad_shape=getDataPointShape(); |
780 |
grad_shape.push_back(functionspace.getDim()); |
781 |
Data out(0.0,grad_shape,functionspace,true); |
782 |
getDomain()->setToGradient(out,*this); |
783 |
blocktimer_increment("grad()", blocktimer_start); |
784 |
return out; |
785 |
} |
786 |
|
787 |
Data |
788 |
Data::grad() const |
789 |
{ |
790 |
if (isEmpty()) |
791 |
{ |
792 |
throw DataException("Error - operation not permitted on instances of DataEmpty."); |
793 |
} |
794 |
return gradOn(escript::function(*getDomain())); |
795 |
} |
796 |
|
797 |
int |
798 |
Data::getDataPointSize() const |
799 |
{ |
800 |
return m_data->getNoValues(); |
801 |
} |
802 |
|
803 |
DataTypes::ValueType::size_type |
804 |
Data::getLength() const |
805 |
{ |
806 |
return m_data->getLength(); |
807 |
} |
808 |
|
809 |
const |
810 |
boost::python::numeric::array |
811 |
Data:: getValueOfDataPoint(int dataPointNo) |
812 |
{ |
813 |
int i, j, k, l; |
814 |
|
815 |
FORCERESOLVE; |
816 |
|
817 |
// |
818 |
// determine the rank and shape of each data point |
819 |
int dataPointRank = getDataPointRank(); |
820 |
const DataTypes::ShapeType& dataPointShape = getDataPointShape(); |
821 |
|
822 |
// |
823 |
// create the numeric array to be returned |
824 |
boost::python::numeric::array numArray(0.0); |
825 |
|
826 |
// |
827 |
// the shape of the returned numeric array will be the same |
828 |
// as that of the data point |
829 |
int arrayRank = dataPointRank; |
830 |
const DataTypes::ShapeType& arrayShape = dataPointShape; |
831 |
|
832 |
// |
833 |
// resize the numeric array to the shape just calculated |
834 |
if (arrayRank==0) { |
835 |
numArray.resize(1); |
836 |
} |
837 |
if (arrayRank==1) { |
838 |
numArray.resize(arrayShape[0]); |
839 |
} |
840 |
if (arrayRank==2) { |
841 |
numArray.resize(arrayShape[0],arrayShape[1]); |
842 |
} |
843 |
if (arrayRank==3) { |
844 |
numArray.resize(arrayShape[0],arrayShape[1],arrayShape[2]); |
845 |
} |
846 |
if (arrayRank==4) { |
847 |
numArray.resize(arrayShape[0],arrayShape[1],arrayShape[2],arrayShape[3]); |
848 |
} |
849 |
|
850 |
if (getNumDataPointsPerSample()>0) { |
851 |
int sampleNo = dataPointNo/getNumDataPointsPerSample(); |
852 |
int dataPointNoInSample = dataPointNo - sampleNo * getNumDataPointsPerSample(); |
853 |
// |
854 |
// Check a valid sample number has been supplied |
855 |
if ((sampleNo >= getNumSamples()) || (sampleNo < 0 )) { |
856 |
throw DataException("Error - Data::convertToNumArray: invalid sampleNo."); |
857 |
} |
858 |
|
859 |
// |
860 |
// Check a valid data point number has been supplied |
861 |
if ((dataPointNoInSample >= getNumDataPointsPerSample()) || (dataPointNoInSample < 0)) { |
862 |
throw DataException("Error - Data::convertToNumArray: invalid dataPointNoInSample."); |
863 |
} |
864 |
// TODO: global error handling |
865 |
// create a view of the data if it is stored locally |
866 |
// DataArrayView dataPointView = getDataPoint(sampleNo, dataPointNoInSample); |
867 |
DataTypes::ValueType::size_type offset=getDataOffset(sampleNo, dataPointNoInSample); |
868 |
|
869 |
|
870 |
switch( dataPointRank ){ |
871 |
case 0 : |
872 |
numArray[0] = getDataAtOffset(offset); |
873 |
break; |
874 |
case 1 : |
875 |
for( i=0; i<dataPointShape[0]; i++ ) |
876 |
numArray[i]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i)); |
877 |
break; |
878 |
case 2 : |
879 |
for( i=0; i<dataPointShape[0]; i++ ) |
880 |
for( j=0; j<dataPointShape[1]; j++) |
881 |
numArray[make_tuple(i,j)]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i,j)); |
882 |
break; |
883 |
case 3 : |
884 |
for( i=0; i<dataPointShape[0]; i++ ) |
885 |
for( j=0; j<dataPointShape[1]; j++ ) |
886 |
for( k=0; k<dataPointShape[2]; k++) |
887 |
numArray[make_tuple(i,j,k)]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i,j,k)); |
888 |
break; |
889 |
case 4 : |
890 |
for( i=0; i<dataPointShape[0]; i++ ) |
891 |
for( j=0; j<dataPointShape[1]; j++ ) |
892 |
for( k=0; k<dataPointShape[2]; k++ ) |
893 |
for( l=0; l<dataPointShape[3]; l++) |
894 |
numArray[make_tuple(i,j,k,l)]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i,j,k,l)); |
895 |
break; |
896 |
} |
897 |
} |
898 |
// |
899 |
// return the array |
900 |
return numArray; |
901 |
|
902 |
} |
903 |
|
904 |
void |
905 |
Data::setValueOfDataPointToPyObject(int dataPointNo, const boost::python::object& py_object) |
906 |
{ |
907 |
// this will throw if the value cannot be represented |
908 |
boost::python::numeric::array num_array(py_object); |
909 |
setValueOfDataPointToArray(dataPointNo,num_array); |
910 |
} |
911 |
|
912 |
void |
913 |
Data::setValueOfDataPointToArray(int dataPointNo, const boost::python::numeric::array& num_array) |
914 |
{ |
915 |
if (isProtected()) { |
916 |
throw DataException("Error - attempt to update protected Data object."); |
917 |
} |
918 |
FORCERESOLVE; |
919 |
// |
920 |
// check rank |
921 |
if (static_cast<unsigned int>(num_array.getrank())<getDataPointRank()) |
922 |
throw DataException("Rank of numarray does not match Data object rank"); |
923 |
|
924 |
// |
925 |
// check shape of num_array |
926 |
for (unsigned int i=0; i<getDataPointRank(); i++) { |
927 |
if (extract<int>(num_array.getshape()[i])!=getDataPointShape()[i]) |
928 |
throw DataException("Shape of numarray does not match Data object rank"); |
929 |
} |
930 |
// |
931 |
// make sure data is expanded: |
932 |
// |
933 |
if (!isExpanded()) { |
934 |
expand(); |
935 |
} |
936 |
if (getNumDataPointsPerSample()>0) { |
937 |
int sampleNo = dataPointNo/getNumDataPointsPerSample(); |
938 |
int dataPointNoInSample = dataPointNo - sampleNo * getNumDataPointsPerSample(); |
939 |
m_data->copyToDataPoint(sampleNo, dataPointNoInSample,num_array); |
940 |
} else { |
941 |
m_data->copyToDataPoint(-1, 0,num_array); |
942 |
} |
943 |
} |
944 |
|
945 |
void |
946 |
Data::setValueOfDataPoint(int dataPointNo, const double value) |
947 |
{ |
948 |
if (isProtected()) { |
949 |
throw DataException("Error - attempt to update protected Data object."); |
950 |
} |
951 |
// |
952 |
// make sure data is expanded: |
953 |
FORCERESOLVE; |
954 |
if (!isExpanded()) { |
955 |
expand(); |
956 |
} |
957 |
if (getNumDataPointsPerSample()>0) { |
958 |
int sampleNo = dataPointNo/getNumDataPointsPerSample(); |
959 |
int dataPointNoInSample = dataPointNo - sampleNo * getNumDataPointsPerSample(); |
960 |
m_data->copyToDataPoint(sampleNo, dataPointNoInSample,value); |
961 |
} else { |
962 |
m_data->copyToDataPoint(-1, 0,value); |
963 |
} |
964 |
} |
965 |
|
966 |
const |
967 |
boost::python::numeric::array |
968 |
Data::getValueOfGlobalDataPoint(int procNo, int dataPointNo) |
969 |
{ |
970 |
size_t length=0; |
971 |
int i, j, k, l, pos; |
972 |
FORCERESOLVE; |
973 |
// |
974 |
// determine the rank and shape of each data point |
975 |
int dataPointRank = getDataPointRank(); |
976 |
const DataTypes::ShapeType& dataPointShape = getDataPointShape(); |
977 |
|
978 |
// |
979 |
// create the numeric array to be returned |
980 |
boost::python::numeric::array numArray(0.0); |
981 |
|
982 |
// |
983 |
// the shape of the returned numeric array will be the same |
984 |
// as that of the data point |
985 |
int arrayRank = dataPointRank; |
986 |
const DataTypes::ShapeType& arrayShape = dataPointShape; |
987 |
|
988 |
// |
989 |
// resize the numeric array to the shape just calculated |
990 |
if (arrayRank==0) { |
991 |
numArray.resize(1); |
992 |
} |
993 |
if (arrayRank==1) { |
994 |
numArray.resize(arrayShape[0]); |
995 |
} |
996 |
if (arrayRank==2) { |
997 |
numArray.resize(arrayShape[0],arrayShape[1]); |
998 |
} |
999 |
if (arrayRank==3) { |
1000 |
numArray.resize(arrayShape[0],arrayShape[1],arrayShape[2]); |
1001 |
} |
1002 |
if (arrayRank==4) { |
1003 |
numArray.resize(arrayShape[0],arrayShape[1],arrayShape[2],arrayShape[3]); |
1004 |
} |
1005 |
|
1006 |
// added for the MPI communication |
1007 |
length=1; |
1008 |
for( i=0; i<arrayRank; i++ ) length *= arrayShape[i]; |
1009 |
double *tmpData = new double[length]; |
1010 |
|
1011 |
// |
1012 |
// load the values for the data point into the numeric array. |
1013 |
|
1014 |
// updated for the MPI case |
1015 |
if( get_MPIRank()==procNo ){ |
1016 |
if (getNumDataPointsPerSample()>0) { |
1017 |
int sampleNo = dataPointNo/getNumDataPointsPerSample(); |
1018 |
int dataPointNoInSample = dataPointNo - sampleNo * getNumDataPointsPerSample(); |
1019 |
// |
1020 |
// Check a valid sample number has been supplied |
1021 |
if ((sampleNo >= getNumSamples()) || (sampleNo < 0 )) { |
1022 |
throw DataException("Error - Data::convertToNumArray: invalid sampleNo."); |
1023 |
} |
1024 |
|
1025 |
// |
1026 |
// Check a valid data point number has been supplied |
1027 |
if ((dataPointNoInSample >= getNumDataPointsPerSample()) || (dataPointNoInSample < 0)) { |
1028 |
throw DataException("Error - Data::convertToNumArray: invalid dataPointNoInSample."); |
1029 |
} |
1030 |
// TODO: global error handling |
1031 |
// create a view of the data if it is stored locally |
1032 |
//DataArrayView dataPointView = getDataPoint(sampleNo, dataPointNoInSample); |
1033 |
DataTypes::ValueType::size_type offset=getDataOffset(sampleNo, dataPointNoInSample); |
1034 |
|
1035 |
// pack the data from the view into tmpData for MPI communication |
1036 |
pos=0; |
1037 |
switch( dataPointRank ){ |
1038 |
case 0 : |
1039 |
tmpData[0] = getDataAtOffset(offset); |
1040 |
break; |
1041 |
case 1 : |
1042 |
for( i=0; i<dataPointShape[0]; i++ ) |
1043 |
tmpData[i]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i)); |
1044 |
break; |
1045 |
case 2 : |
1046 |
for( i=0; i<dataPointShape[0]; i++ ) |
1047 |
for( j=0; j<dataPointShape[1]; j++, pos++ ) |
1048 |
tmpData[pos]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i,j)); |
1049 |
break; |
1050 |
case 3 : |
1051 |
for( i=0; i<dataPointShape[0]; i++ ) |
1052 |
for( j=0; j<dataPointShape[1]; j++ ) |
1053 |
for( k=0; k<dataPointShape[2]; k++, pos++ ) |
1054 |
tmpData[pos]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i,j,k)); |
1055 |
break; |
1056 |
case 4 : |
1057 |
for( i=0; i<dataPointShape[0]; i++ ) |
1058 |
for( j=0; j<dataPointShape[1]; j++ ) |
1059 |
for( k=0; k<dataPointShape[2]; k++ ) |
1060 |
for( l=0; l<dataPointShape[3]; l++, pos++ ) |
1061 |
tmpData[pos]=getDataAtOffset(offset+DataTypes::getRelIndex(dataPointShape, i,j,k,l)); |
1062 |
break; |
1063 |
} |
1064 |
} |
1065 |
} |
1066 |
#ifdef PASO_MPI |
1067 |
// broadcast the data to all other processes |
1068 |
MPI_Bcast( tmpData, length, MPI_DOUBLE, procNo, get_MPIComm() ); |
1069 |
#endif |
1070 |
|
1071 |
// unpack the data |
1072 |
switch( dataPointRank ){ |
1073 |
case 0 : |
1074 |
numArray[0]=tmpData[0]; |
1075 |
break; |
1076 |
case 1 : |
1077 |
for( i=0; i<dataPointShape[0]; i++ ) |
1078 |
numArray[i]=tmpData[i]; |
1079 |
break; |
1080 |
case 2 : |
1081 |
for( i=0; i<dataPointShape[0]; i++ ) |
1082 |
for( j=0; j<dataPointShape[1]; j++ ) |
1083 |
numArray[make_tuple(i,j)]=tmpData[i+j*dataPointShape[0]]; |
1084 |
break; |
1085 |
case 3 : |
1086 |
for( i=0; i<dataPointShape[0]; i++ ) |
1087 |
for( j=0; j<dataPointShape[1]; j++ ) |
1088 |
for( k=0; k<dataPointShape[2]; k++ ) |
1089 |
numArray[make_tuple(i,j,k)]=tmpData[i+dataPointShape[0]*(j*+k*dataPointShape[1])]; |
1090 |
break; |
1091 |
case 4 : |
1092 |
for( i=0; i<dataPointShape[0]; i++ ) |
1093 |
for( j=0; j<dataPointShape[1]; j++ ) |
1094 |
for( k=0; k<dataPointShape[2]; k++ ) |
1095 |
for( l=0; l<dataPointShape[3]; l++ ) |
1096 |
numArray[make_tuple(i,j,k,l)]=tmpData[i+dataPointShape[0]*(j*+dataPointShape[1]*(k+l*dataPointShape[2]))]; |
1097 |
break; |
1098 |
} |
1099 |
|
1100 |
delete [] tmpData; |
1101 |
// |
1102 |
// return the loaded array |
1103 |
return numArray; |
1104 |
} |
1105 |
|
1106 |
|
1107 |
boost::python::numeric::array |
1108 |
Data::integrate_const() const |
1109 |
{ |
1110 |
if (isLazy()) |
1111 |
{ |
1112 |
throw DataException("Error - cannot integrate for constant lazy data."); |
1113 |
} |
1114 |
return integrateWorker(); |
1115 |
} |
1116 |
|
1117 |
boost::python::numeric::array |
1118 |
Data::integrate() |
1119 |
{ |
1120 |
if (isLazy()) |
1121 |
{ |
1122 |
expand(); |
1123 |
} |
1124 |
return integrateWorker(); |
1125 |
} |
1126 |
|
1127 |
|
1128 |
|
1129 |
boost::python::numeric::array |
1130 |
Data::integrateWorker() const |
1131 |
{ |
1132 |
int index; |
1133 |
int rank = getDataPointRank(); |
1134 |
DataTypes::ShapeType shape = getDataPointShape(); |
1135 |
int dataPointSize = getDataPointSize(); |
1136 |
|
1137 |
// |
1138 |
// calculate the integral values |
1139 |
vector<double> integrals(dataPointSize); |
1140 |
vector<double> integrals_local(dataPointSize); |
1141 |
const AbstractContinuousDomain* dom=dynamic_cast<const AbstractContinuousDomain*>(getDomain().get()); |
1142 |
if (dom==0) |
1143 |
{ |
1144 |
throw DataException("Can not integrate over non-continuous domains."); |
1145 |
} |
1146 |
#ifdef PASO_MPI |
1147 |
dom->setToIntegrals(integrals_local,*this); |
1148 |
// Global sum: use an array instead of a vector because elements of array are guaranteed to be contiguous in memory |
1149 |
double *tmp = new double[dataPointSize]; |
1150 |
double *tmp_local = new double[dataPointSize]; |
1151 |
for (int i=0; i<dataPointSize; i++) { tmp_local[i] = integrals_local[i]; } |
1152 |
MPI_Allreduce( &tmp_local[0], &tmp[0], dataPointSize, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD ); |
1153 |
for (int i=0; i<dataPointSize; i++) { integrals[i] = tmp[i]; } |
1154 |
delete[] tmp; |
1155 |
delete[] tmp_local; |
1156 |
#else |
1157 |
dom->setToIntegrals(integrals,*this); |
1158 |
#endif |
1159 |
|
1160 |
// |
1161 |
// create the numeric array to be returned |
1162 |
// and load the array with the integral values |
1163 |
boost::python::numeric::array bp_array(1.0); |
1164 |
if (rank==0) { |
1165 |
bp_array.resize(1); |
1166 |
index = 0; |
1167 |
bp_array[0] = integrals[index]; |
1168 |
} |
1169 |
if (rank==1) { |
1170 |
bp_array.resize(shape[0]); |
1171 |
for (int i=0; i<shape[0]; i++) { |
1172 |
index = i; |
1173 |
bp_array[i] = integrals[index]; |
1174 |
} |
1175 |
} |
1176 |
if (rank==2) { |
1177 |
bp_array.resize(shape[0],shape[1]); |
1178 |
for (int i=0; i<shape[0]; i++) { |
1179 |
for (int j=0; j<shape[1]; j++) { |
1180 |
index = i + shape[0] * j; |
1181 |
bp_array[make_tuple(i,j)] = integrals[index]; |
1182 |
} |
1183 |
} |
1184 |
} |
1185 |
if (rank==3) { |
1186 |
bp_array.resize(shape[0],shape[1],shape[2]); |
1187 |
for (int i=0; i<shape[0]; i++) { |
1188 |
for (int j=0; j<shape[1]; j++) { |
1189 |
for (int k=0; k<shape[2]; k++) { |
1190 |
index = i + shape[0] * ( j + shape[1] * k ); |
1191 |
bp_array[make_tuple(i,j,k)] = integrals[index]; |
1192 |
} |
1193 |
} |
1194 |
} |
1195 |
} |
1196 |
if (rank==4) { |
1197 |
bp_array.resize(shape[0],shape[1],shape[2],shape[3]); |
1198 |
for (int i=0; i<shape[0]; i++) { |
1199 |
for (int j=0; j<shape[1]; j++) { |
1200 |
for (int k=0; k<shape[2]; k++) { |
1201 |
for (int l=0; l<shape[3]; l++) { |
1202 |
index = i + shape[0] * ( j + shape[1] * ( k + shape[2] * l ) ); |
1203 |
bp_array[make_tuple(i,j,k,l)] = integrals[index]; |
1204 |
} |
1205 |
} |
1206 |
} |
1207 |
} |
1208 |
} |
1209 |
|
1210 |
// |
1211 |
// return the loaded array |
1212 |
return bp_array; |
1213 |
} |
1214 |
|
1215 |
Data |
1216 |
Data::sin() const |
1217 |
{ |
1218 |
if (isLazy()) |
1219 |
{ |
1220 |
DataLazy* c=new DataLazy(borrowDataPtr(),SIN); |
1221 |
return Data(c); |
1222 |
} |
1223 |
return C_TensorUnaryOperation<double (*)(double)>(*this, ::sin); |
1224 |
} |
1225 |
|
1226 |
Data |
1227 |
Data::cos() const |
1228 |
{ |
1229 |
if (isLazy()) |
1230 |
{ |
1231 |
DataLazy* c=new DataLazy(borrowDataPtr(),COS); |
1232 |
return Data(c); |
1233 |
} |
1234 |
return C_TensorUnaryOperation<double (*)(double)>(*this, ::cos); |
1235 |
} |
1236 |
|
1237 |
Data |
1238 |
Data::tan() const |
1239 |
{ |
1240 |
if (isLazy()) |
1241 |
{ |
1242 |
DataLazy* c=new DataLazy(borrowDataPtr(),TAN); |
1243 |
return Data(c); |
1244 |
} |
1245 |
return C_TensorUnaryOperation<double (*)(double)>(*this, ::tan); |
1246 |
} |
1247 |
|
1248 |
Data |
1249 |
Data::asin() const |
1250 |
{ |
1251 |
if (isLazy()) |
1252 |
{ |
1253 |
DataLazy* c=new DataLazy(borrowDataPtr(),ASIN); |
1254 |
return Data(c); |
1255 |
} |
1256 |
return C_TensorUnaryOperation<double (*)(double)>(*this, ::asin); |
1257 |
} |
1258 |
|
1259 |
Data |
1260 |
Data::acos() const |
1261 |
{ |
1262 |
if (isLazy()) |
1263 |
{ |
1264 |
DataLazy* c=new DataLazy(borrowDataPtr(),ACOS); |
1265 |
return Data(c); |
1266 |
} |
1267 |
return C_TensorUnaryOperation<double (*)(double)>(*this, ::acos); |
1268 |
} |
1269 |
|
1270 |
|
1271 |
Data |
1272 |
Data::atan() const |
1273 |
{ |
1274 |
if (isLazy()) |
1275 |
{ |
1276 |
DataLazy* c=new DataLazy(borrowDataPtr(),ATAN); |
1277 |
return Data(c); |
1278 |
} |
1279 |
return C_TensorUnaryOperation<double (*)(double)>(*this, ::atan); |
1280 |
} |
1281 |
|
1282 |
Data |
1283 |
Data::sinh() const |
1284 |
{ |
1285 |
if (isLazy()) |
1286 |
{ |
1287 |
DataLazy* c=new DataLazy(borrowDataPtr(),SINH); |
1288 |
return Data(c); |
1289 |
} |
1290 |
return C_TensorUnaryOperation<double (*)(double)>(*this, ::sinh); |
1291 |
} |
1292 |
|
1293 |
Data |
1294 |
Data::cosh() const |
1295 |
{ |
1296 |
if (isLazy()) |
1297 |
{ |
1298 |
DataLazy* c=new DataLazy(borrowDataPtr(),COSH); |
1299 |
return Data(c); |
1300 |
} |
1301 |
return C_TensorUnaryOperation<double (*)(double)>(*this, ::cosh); |
1302 |
} |
1303 |
|
1304 |
Data |
1305 |
Data::tanh() const |
1306 |
{ |
1307 |
if (isLazy()) |
1308 |
{ |
1309 |
DataLazy* c=new DataLazy(borrowDataPtr(),TANH); |
1310 |
return Data(c); |
1311 |
} |
1312 |
return C_TensorUnaryOperation<double (*)(double)>(*this, ::tanh); |
1313 |
} |
1314 |
|
1315 |
|
1316 |
Data |
1317 |
Data::erf() const |
1318 |
{ |
1319 |
#if defined (_WIN32) && !defined(__INTEL_COMPILER) |
1320 |
throw DataException("Error - Data:: erf function is not supported on _WIN32 platforms."); |
1321 |
#else |
1322 |
if (isLazy()) |
1323 |
{ |
1324 |
DataLazy* c=new DataLazy(borrowDataPtr(),ERF); |
1325 |
return Data(c); |
1326 |
} |
1327 |
return C_TensorUnaryOperation(*this, ::erf); |
1328 |
#endif |
1329 |
} |
1330 |
|
1331 |
Data |
1332 |
Data::asinh() const |
1333 |
{ |
1334 |
if (isLazy()) |
1335 |
{ |
1336 |
DataLazy* c=new DataLazy(borrowDataPtr(),ASINH); |
1337 |
return Data(c); |
1338 |
} |
1339 |
#if defined (_WIN32) && !defined(__INTEL_COMPILER) |
1340 |
return C_TensorUnaryOperation(*this, escript::asinh_substitute); |
1341 |
#else |
1342 |
return C_TensorUnaryOperation(*this, ::asinh); |
1343 |
#endif |
1344 |
} |
1345 |
|
1346 |
Data |
1347 |
Data::acosh() const |
1348 |
{ |
1349 |
if (isLazy()) |
1350 |
{ |
1351 |
DataLazy* c=new DataLazy(borrowDataPtr(),ACOSH); |
1352 |
return Data(c); |
1353 |
} |
1354 |
#if defined (_WIN32) && !defined(__INTEL_COMPILER) |
1355 |
return C_TensorUnaryOperation(*this, escript::acosh_substitute); |
1356 |
#else |
1357 |
return C_TensorUnaryOperation(*this, ::acosh); |
1358 |
#endif |
1359 |
} |
1360 |
|
1361 |
Data |
1362 |
Data::atanh() const |
1363 |
{ |
1364 |
if (isLazy()) |
1365 |
{ |
1366 |
DataLazy* c=new DataLazy(borrowDataPtr(),ATANH); |
1367 |
return Data(c); |
1368 |
} |
1369 |
#if defined (_WIN32) && !defined(__INTEL_COMPILER) |
1370 |
return C_TensorUnaryOperation(*this, escript::atanh_substitute); |
1371 |
#else |
1372 |
return C_TensorUnaryOperation(*this, ::atanh); |
1373 |
#endif |
1374 |
} |
1375 |
|
1376 |
Data |
1377 |
Data::log10() const |
1378 |
{ if (isLazy()) |
1379 |
{ |
1380 |
DataLazy* c=new DataLazy(borrowDataPtr(),LOG10); |
1381 |
return Data(c); |
1382 |
} |
1383 |
return C_TensorUnaryOperation<double (*)(double)>(*this, ::log10); |
1384 |
} |
1385 |
|
1386 |
Data |
1387 |
Data::log() const |
1388 |
{ |
1389 |
if (isLazy()) |
1390 |
{ |
1391 |
DataLazy* c=new DataLazy(borrowDataPtr(),LOG); |
1392 |
return Data(c); |
1393 |
} |
1394 |
return C_TensorUnaryOperation<double (*)(double)>(*this, ::log); |
1395 |
} |
1396 |
|
1397 |
Data |
1398 |
Data::sign() const |
1399 |
{ |
1400 |
if (isLazy()) |
1401 |
{ |
1402 |
DataLazy* c=new DataLazy(borrowDataPtr(),SIGN); |
1403 |
return Data(c); |
1404 |
} |
1405 |
return C_TensorUnaryOperation(*this, escript::fsign); |
1406 |
} |
1407 |
|
1408 |
Data |
1409 |
Data::abs() const |
1410 |
{ |
1411 |
if (isLazy()) |
1412 |
{ |
1413 |
DataLazy* c=new DataLazy(borrowDataPtr(),ABS); |
1414 |
return Data(c); |
1415 |
} |
1416 |
return C_TensorUnaryOperation<double (*)(double)>(*this, ::fabs); |
1417 |
} |
1418 |
|
1419 |
Data |
1420 |
Data::neg() const |
1421 |
{ |
1422 |
if (isLazy()) |
1423 |
{ |
1424 |
DataLazy* c=new DataLazy(borrowDataPtr(),NEG); |
1425 |
return Data(c); |
1426 |
} |
1427 |
return C_TensorUnaryOperation(*this, negate<double>()); |
1428 |
} |
1429 |
|
1430 |
Data |
1431 |
Data::pos() const |
1432 |
{ |
1433 |
// not doing lazy check here is deliberate. |
1434 |
// since a deep copy of lazy data should be cheap, I'll just let it happen now |
1435 |
Data result; |
1436 |
// perform a deep copy |
1437 |
result.copy(*this); |
1438 |
return result; |
1439 |
} |
1440 |
|
1441 |
Data |
1442 |
Data::exp() const |
1443 |
{ |
1444 |
if (isLazy()) |
1445 |
{ |
1446 |
DataLazy* c=new DataLazy(borrowDataPtr(),EXP); |
1447 |
return Data(c); |
1448 |
} |
1449 |
return C_TensorUnaryOperation<double (*)(double)>(*this, ::exp); |
1450 |
} |
1451 |
|
1452 |
Data |
1453 |
Data::sqrt() const |
1454 |
{ |
1455 |
if (isLazy()) |
1456 |
{ |
1457 |
DataLazy* c=new DataLazy(borrowDataPtr(),SQRT); |
1458 |
return Data(c); |
1459 |
} |
1460 |
return C_TensorUnaryOperation<double (*)(double)>(*this, ::sqrt); |
1461 |
} |
1462 |
|
1463 |
double |
1464 |
Data::Lsup_const() const |
1465 |
{ |
1466 |
if (isLazy()) |
1467 |
{ |
1468 |
throw DataException("Error - cannot compute Lsup for constant lazy data."); |
1469 |
} |
1470 |
return LsupWorker(); |
1471 |
} |
1472 |
|
1473 |
double |
1474 |
Data::Lsup() |
1475 |
{ |
1476 |
if (isLazy()) |
1477 |
{ |
1478 |
resolve(); |
1479 |
} |
1480 |
return LsupWorker(); |
1481 |
} |
1482 |
|
1483 |
double |
1484 |
Data::sup_const() const |
1485 |
{ |
1486 |
if (isLazy()) |
1487 |
{ |
1488 |
throw DataException("Error - cannot compute sup for constant lazy data."); |
1489 |
} |
1490 |
return supWorker(); |
1491 |
} |
1492 |
|
1493 |
double |
1494 |
Data::sup() |
1495 |
{ |
1496 |
if (isLazy()) |
1497 |
{ |
1498 |
resolve(); |
1499 |
} |
1500 |
return supWorker(); |
1501 |
} |
1502 |
|
1503 |
double |
1504 |
Data::inf_const() const |
1505 |
{ |
1506 |
if (isLazy()) |
1507 |
{ |
1508 |
throw DataException("Error - cannot compute inf for constant lazy data."); |
1509 |
} |
1510 |
return infWorker(); |
1511 |
} |
1512 |
|
1513 |
double |
1514 |
Data::inf() |
1515 |
{ |
1516 |
if (isLazy()) |
1517 |
{ |
1518 |
resolve(); |
1519 |
} |
1520 |
return infWorker(); |
1521 |
} |
1522 |
|
1523 |
double |
1524 |
Data::LsupWorker() const |
1525 |
{ |
1526 |
double localValue; |
1527 |
// |
1528 |
// set the initial absolute maximum value to zero |
1529 |
|
1530 |
AbsMax abs_max_func; |
1531 |
localValue = algorithm(abs_max_func,0); |
1532 |
#ifdef PASO_MPI |
1533 |
double globalValue; |
1534 |
MPI_Allreduce( &localValue, &globalValue, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD ); |
1535 |
return globalValue; |
1536 |
#else |
1537 |
return localValue; |
1538 |
#endif |
1539 |
} |
1540 |
|
1541 |
double |
1542 |
Data::supWorker() const |
1543 |
{ |
1544 |
double localValue; |
1545 |
// |
1546 |
// set the initial maximum value to min possible double |
1547 |
FMax fmax_func; |
1548 |
localValue = algorithm(fmax_func,numeric_limits<double>::max()*-1); |
1549 |
#ifdef PASO_MPI |
1550 |
double globalValue; |
1551 |
MPI_Allreduce( &localValue, &globalValue, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD ); |
1552 |
return globalValue; |
1553 |
#else |
1554 |
return localValue; |
1555 |
#endif |
1556 |
} |
1557 |
|
1558 |
double |
1559 |
Data::infWorker() const |
1560 |
{ |
1561 |
double localValue; |
1562 |
// |
1563 |
// set the initial minimum value to max possible double |
1564 |
FMin fmin_func; |
1565 |
localValue = algorithm(fmin_func,numeric_limits<double>::max()); |
1566 |
#ifdef PASO_MPI |
1567 |
double globalValue; |
1568 |
MPI_Allreduce( &localValue, &globalValue, 1, MPI_DOUBLE, MPI_MIN, MPI_COMM_WORLD ); |
1569 |
return globalValue; |
1570 |
#else |
1571 |
return localValue; |
1572 |
#endif |
1573 |
} |
1574 |
|
1575 |
/* TODO */ |
1576 |
/* global reduction */ |
1577 |
Data |
1578 |
Data::maxval() const |
1579 |
{ |
1580 |
if (isLazy()) |
1581 |
{ |
1582 |
Data temp(*this); // to get around the fact that you can't resolve a const Data |
1583 |
temp.resolve(); |
1584 |
return temp.maxval(); |
1585 |
} |
1586 |
// |
1587 |
// set the initial maximum value to min possible double |
1588 |
FMax fmax_func; |
1589 |
return dp_algorithm(fmax_func,numeric_limits<double>::max()*-1); |
1590 |
} |
1591 |
|
1592 |
Data |
1593 |
Data::minval() const |
1594 |
{ |
1595 |
if (isLazy()) |
1596 |
{ |
1597 |
Data temp(*this); // to get around the fact that you can't resolve a const Data |
1598 |
temp.resolve(); |
1599 |
return temp.minval(); |
1600 |
} |
1601 |
// |
1602 |
// set the initial minimum value to max possible double |
1603 |
FMin fmin_func; |
1604 |
return dp_algorithm(fmin_func,numeric_limits<double>::max()); |
1605 |
} |
1606 |
|
1607 |
Data |
1608 |
Data::swapaxes(const int axis0, const int axis1) const |
1609 |
{ |
1610 |
if (isLazy()) |
1611 |
{ |
1612 |
Data temp(*this); |
1613 |
temp.resolve(); |
1614 |
return temp.swapaxes(axis0,axis1); |
1615 |
} |
1616 |
int axis0_tmp,axis1_tmp; |
1617 |
DataTypes::ShapeType s=getDataPointShape(); |
1618 |
DataTypes::ShapeType ev_shape; |
1619 |
// Here's the equivalent of python s_out=s[axis_offset:]+s[:axis_offset] |
1620 |
// which goes thru all shape vector elements starting with axis_offset (at index=rank wrap around to 0) |
1621 |
int rank=getDataPointRank(); |
1622 |
if (rank<2) { |
1623 |
throw DataException("Error - Data::swapaxes argument must have at least rank 2."); |
1624 |
} |
1625 |
if (axis0<0 || axis0>rank-1) { |
1626 |
throw DataException("Error - Data::swapaxes: axis0 must be between 0 and rank-1=" + rank-1); |
1627 |
} |
1628 |
if (axis1<0 || axis1>rank-1) { |
1629 |
throw DataException("Error - Data::swapaxes: axis1 must be between 0 and rank-1=" + rank-1); |
1630 |
} |
1631 |
if (axis0 == axis1) { |
1632 |
throw DataException("Error - Data::swapaxes: axis indices must be different."); |
1633 |
} |
1634 |
if (axis0 > axis1) { |
1635 |
axis0_tmp=axis1; |
1636 |
axis1_tmp=axis0; |
1637 |
} else { |
1638 |
axis0_tmp=axis0; |
1639 |
axis1_tmp=axis1; |
1640 |
} |
1641 |
for (int i=0; i<rank; i++) { |
1642 |
if (i == axis0_tmp) { |
1643 |
ev_shape.push_back(s[axis1_tmp]); |
1644 |
} else if (i == axis1_tmp) { |
1645 |
ev_shape.push_back(s[axis0_tmp]); |
1646 |
} else { |
1647 |
ev_shape.push_back(s[i]); |
1648 |
} |
1649 |
} |
1650 |
Data ev(0.,ev_shape,getFunctionSpace()); |
1651 |
ev.typeMatchRight(*this); |
1652 |
m_data->swapaxes(ev.m_data.get(), axis0_tmp, axis1_tmp); |
1653 |
return ev; |
1654 |
|
1655 |
} |
1656 |
|
1657 |
Data |
1658 |
Data::symmetric() const |
1659 |
{ |
1660 |
// check input |
1661 |
DataTypes::ShapeType s=getDataPointShape(); |
1662 |
if (getDataPointRank()==2) { |
1663 |
if(s[0] != s[1]) |
1664 |
throw DataException("Error - Data::symmetric can only be calculated for rank 2 object with equal first and second dimension."); |
1665 |
} |
1666 |
else if (getDataPointRank()==4) { |
1667 |
if(!(s[0] == s[2] && s[1] == s[3])) |
1668 |
throw DataException("Error - Data::symmetric can only be calculated for rank 4 object with dim0==dim2 and dim1==dim3."); |
1669 |
} |
1670 |
else { |
1671 |
throw DataException("Error - Data::symmetric can only be calculated for rank 2 or 4 object."); |
1672 |
} |
1673 |
if (isLazy()) |
1674 |
{ |
1675 |
DataLazy* c=new DataLazy(borrowDataPtr(),SYM); |
1676 |
return Data(c); |
1677 |
} |
1678 |
Data ev(0.,getDataPointShape(),getFunctionSpace()); |
1679 |
ev.typeMatchRight(*this); |
1680 |
m_data->symmetric(ev.m_data.get()); |
1681 |
return ev; |
1682 |
} |
1683 |
|
1684 |
Data |
1685 |
Data::nonsymmetric() const |
1686 |
{ |
1687 |
if (isLazy()) |
1688 |
{ |
1689 |
DataLazy* c=new DataLazy(borrowDataPtr(),NSYM); |
1690 |
return Data(c); |
1691 |
} |
1692 |
// check input |
1693 |
DataTypes::ShapeType s=getDataPointShape(); |
1694 |
if (getDataPointRank()==2) { |
1695 |
if(s[0] != s[1]) |
1696 |
throw DataException("Error - Data::nonsymmetric can only be calculated for rank 2 object with equal first and second dimension."); |
1697 |
DataTypes::ShapeType ev_shape; |
1698 |
ev_shape.push_back(s[0]); |
1699 |
ev_shape.push_back(s[1]); |
1700 |
Data ev(0.,ev_shape,getFunctionSpace()); |
1701 |
ev.typeMatchRight(*this); |
1702 |
m_data->nonsymmetric(ev.m_data.get()); |
1703 |
return ev; |
1704 |
} |
1705 |
else if (getDataPointRank()==4) { |
1706 |
if(!(s[0] == s[2] && s[1] == s[3])) |
1707 |
throw DataException("Error - Data::nonsymmetric can only be calculated for rank 4 object with dim0==dim2 and dim1==dim3."); |
1708 |
DataTypes::ShapeType ev_shape; |
1709 |
ev_shape.push_back(s[0]); |
1710 |
ev_shape.push_back(s[1]); |
1711 |
ev_shape.push_back(s[2]); |
1712 |
ev_shape.push_back(s[3]); |
1713 |
Data ev(0.,ev_shape,getFunctionSpace()); |
1714 |
ev.typeMatchRight(*this); |
1715 |
m_data->nonsymmetric(ev.m_data.get()); |
1716 |
return ev; |
1717 |
} |
1718 |
else { |
1719 |
throw DataException("Error - Data::nonsymmetric can only be calculated for rank 2 or 4 object."); |
1720 |
} |
1721 |
} |
1722 |
|
1723 |
Data |
1724 |
Data::trace(int axis_offset) const |
1725 |
{ |
1726 |
if (isLazy()) |
1727 |
{ |
1728 |
DataLazy* c=new DataLazy(borrowDataPtr(),TRACE,axis_offset); |
1729 |
return Data(c); |
1730 |
} |
1731 |
DataTypes::ShapeType s=getDataPointShape(); |
1732 |
if (getDataPointRank()==2) { |
1733 |
DataTypes::ShapeType ev_shape; |
1734 |
Data ev(0.,ev_shape,getFunctionSpace()); |
1735 |
ev.typeMatchRight(*this); |
1736 |
m_data->trace(ev.m_data.get(), axis_offset); |
1737 |
return ev; |
1738 |
} |
1739 |
if (getDataPointRank()==3) { |
1740 |
DataTypes::ShapeType ev_shape; |
1741 |
if (axis_offset==0) { |
1742 |
int s2=s[2]; |
1743 |
ev_shape.push_back(s2); |
1744 |
} |
1745 |
else if (axis_offset==1) { |
1746 |
int s0=s[0]; |
1747 |
ev_shape.push_back(s0); |
1748 |
} |
1749 |
Data ev(0.,ev_shape,getFunctionSpace()); |
1750 |
ev.typeMatchRight(*this); |
1751 |
m_data->trace(ev.m_data.get(), axis_offset); |
1752 |
return ev; |
1753 |
} |
1754 |
if (getDataPointRank()==4) { |
1755 |
DataTypes::ShapeType ev_shape; |
1756 |
if (axis_offset==0) { |
1757 |
ev_shape.push_back(s[2]); |
1758 |
ev_shape.push_back(s[3]); |
1759 |
} |
1760 |
else if (axis_offset==1) { |
1761 |
ev_shape.push_back(s[0]); |
1762 |
ev_shape.push_back(s[3]); |
1763 |
} |
1764 |
else if (axis_offset==2) { |
1765 |
ev_shape.push_back(s[0]); |
1766 |
ev_shape.push_back(s[1]); |
1767 |
} |
1768 |
Data ev(0.,ev_shape,getFunctionSpace()); |
1769 |
ev.typeMatchRight(*this); |
1770 |
m_data->trace(ev.m_data.get(), axis_offset); |
1771 |
return ev; |
1772 |
} |
1773 |
else { |
1774 |
throw DataException("Error - Data::trace can only be calculated for rank 2, 3 or 4 object."); |
1775 |
} |
1776 |
} |
1777 |
|
1778 |
Data |
1779 |
Data::transpose(int axis_offset) const |
1780 |
{ |
1781 |
if (isLazy()) |
1782 |
{ |
1783 |
DataLazy* c=new DataLazy(borrowDataPtr(),TRANS,axis_offset); |
1784 |
return Data(c); |
1785 |
} |
1786 |
DataTypes::ShapeType s=getDataPointShape(); |
1787 |
DataTypes::ShapeType ev_shape; |
1788 |
// Here's the equivalent of python s_out=s[axis_offset:]+s[:axis_offset] |
1789 |
// which goes thru all shape vector elements starting with axis_offset (at index=rank wrap around to 0) |
1790 |
int rank=getDataPointRank(); |
1791 |
if (axis_offset<0 || axis_offset>rank) { |
1792 |
throw DataException("Error - Data::transpose must have 0 <= axis_offset <= rank=" + rank); |
1793 |
} |
1794 |
for (int i=0; i<rank; i++) { |
1795 |
int index = (axis_offset+i)%rank; |
1796 |
ev_shape.push_back(s[index]); // Append to new shape |
1797 |
} |
1798 |
Data ev(0.,ev_shape,getFunctionSpace()); |
1799 |
ev.typeMatchRight(*this); |
1800 |
m_data->transpose(ev.m_data.get(), axis_offset); |
1801 |
return ev; |
1802 |
} |
1803 |
|
1804 |
Data |
1805 |
Data::eigenvalues() const |
1806 |
{ |
1807 |
if (isLazy()) |
1808 |
{ |
1809 |
Data temp(*this); // to get around the fact that you can't resolve a const Data |
1810 |
temp.resolve(); |
1811 |
return temp.eigenvalues(); |
1812 |
} |
1813 |
// check input |
1814 |
DataTypes::ShapeType s=getDataPointShape(); |
1815 |
if (getDataPointRank()!=2) |
1816 |
throw DataException("Error - Data::eigenvalues can only be calculated for rank 2 object."); |
1817 |
if(s[0] != s[1]) |
1818 |
throw DataException("Error - Data::eigenvalues can only be calculated for object with equal first and second dimension."); |
1819 |
// create return |
1820 |
DataTypes::ShapeType ev_shape(1,s[0]); |
1821 |
Data ev(0.,ev_shape,getFunctionSpace()); |
1822 |
ev.typeMatchRight(*this); |
1823 |
m_data->eigenvalues(ev.m_data.get()); |
1824 |
return ev; |
1825 |
} |
1826 |
|
1827 |
const boost::python::tuple |
1828 |
Data::eigenvalues_and_eigenvectors(const double tol) const |
1829 |
{ |
1830 |
if (isLazy()) |
1831 |
{ |
1832 |
Data temp(*this); // to get around the fact that you can't resolve a const Data |
1833 |
temp.resolve(); |
1834 |
return temp.eigenvalues_and_eigenvectors(tol); |
1835 |
} |
1836 |
DataTypes::ShapeType s=getDataPointShape(); |
1837 |
if (getDataPointRank()!=2) |
1838 |
throw DataException("Error - Data::eigenvalues and eigenvectors can only be calculated for rank 2 object."); |
1839 |
if(s[0] != s[1]) |
1840 |
throw DataException("Error - Data::eigenvalues and eigenvectors can only be calculated for object with equal first and second dimension."); |
1841 |
// create return |
1842 |
DataTypes::ShapeType ev_shape(1,s[0]); |
1843 |
Data ev(0.,ev_shape,getFunctionSpace()); |
1844 |
ev.typeMatchRight(*this); |
1845 |
DataTypes::ShapeType V_shape(2,s[0]); |
1846 |
Data V(0.,V_shape,getFunctionSpace()); |
1847 |
V.typeMatchRight(*this); |
1848 |
m_data->eigenvalues_and_eigenvectors(ev.m_data.get(),V.m_data.get(),tol); |
1849 |
return make_tuple(boost::python::object(ev),boost::python::object(V)); |
1850 |
} |
1851 |
|
1852 |
const boost::python::tuple |
1853 |
Data::minGlobalDataPoint() const |
1854 |
{ |
1855 |
// NB: calc_minGlobalDataPoint( had to be split off from minGlobalDataPoint( as boost::make_tuple causes an |
1856 |
// abort (for unknown reasons) if there are openmp directives with it in the |
1857 |
// surrounding function |
1858 |
|
1859 |
int DataPointNo; |
1860 |
int ProcNo; |
1861 |
calc_minGlobalDataPoint(ProcNo,DataPointNo); |
1862 |
return make_tuple(ProcNo,DataPointNo); |
1863 |
} |
1864 |
|
1865 |
void |
1866 |
Data::calc_minGlobalDataPoint(int& ProcNo, |
1867 |
int& DataPointNo) const |
1868 |
{ |
1869 |
if (isLazy()) |
1870 |
{ |
1871 |
Data temp(*this); // to get around the fact that you can't resolve a const Data |
1872 |
temp.resolve(); |
1873 |
return temp.calc_minGlobalDataPoint(ProcNo,DataPointNo); |
1874 |
} |
1875 |
int i,j; |
1876 |
int lowi=0,lowj=0; |
1877 |
double min=numeric_limits<double>::max(); |
1878 |
|
1879 |
Data temp=minval(); |
1880 |
|
1881 |
int numSamples=temp.getNumSamples(); |
1882 |
int numDPPSample=temp.getNumDataPointsPerSample(); |
1883 |
|
1884 |
double next,local_min; |
1885 |
int local_lowi=0,local_lowj=0; |
1886 |
|
1887 |
#pragma omp parallel firstprivate(local_lowi,local_lowj) private(next,local_min) |
1888 |
{ |
1889 |
local_min=min; |
1890 |
#pragma omp for private(i,j) schedule(static) |
1891 |
for (i=0; i<numSamples; i++) { |
1892 |
for (j=0; j<numDPPSample; j++) { |
1893 |
next=temp.getDataAtOffset(temp.getDataOffset(i,j)); |
1894 |
if (next<local_min) { |
1895 |
local_min=next; |
1896 |
local_lowi=i; |
1897 |
local_lowj=j; |
1898 |
} |
1899 |
} |
1900 |
} |
1901 |
#pragma omp critical |
1902 |
if (local_min<min) { |
1903 |
min=local_min; |
1904 |
lowi=local_lowi; |
1905 |
lowj=local_lowj; |
1906 |
} |
1907 |
} |
1908 |
|
1909 |
#ifdef PASO_MPI |
1910 |
// determine the processor on which the minimum occurs |
1911 |
next = temp.getDataPoint(lowi,lowj); |
1912 |
int lowProc = 0; |
1913 |
double *globalMins = new double[get_MPISize()+1]; |
1914 |
int error; |
1915 |
error = MPI_Gather ( &next, 1, MPI_DOUBLE, globalMins, 1, MPI_DOUBLE, 0, get_MPIComm() ); |
1916 |
|
1917 |
if( get_MPIRank()==0 ){ |
1918 |
next = globalMins[lowProc]; |
1919 |
for( i=1; i<get_MPISize(); i++ ) |
1920 |
if( next>globalMins[i] ){ |
1921 |
lowProc = i; |
1922 |
next = globalMins[i]; |
1923 |
} |
1924 |
} |
1925 |
MPI_Bcast( &lowProc, 1, MPI_DOUBLE, 0, get_MPIComm() ); |
1926 |
|
1927 |
delete [] globalMins; |
1928 |
ProcNo = lowProc; |
1929 |
#else |
1930 |
ProcNo = 0; |
1931 |
#endif |
1932 |
DataPointNo = lowj + lowi * numDPPSample; |
1933 |
} |
1934 |
|
1935 |
void |
1936 |
Data::saveDX(std::string fileName) const |
1937 |
{ |
1938 |
if (isEmpty()) |
1939 |
{ |
1940 |
throw DataException("Error - Operations not permitted on instances of DataEmpty."); |
1941 |
} |
1942 |
if (isLazy()) |
1943 |
{ |
1944 |
Data temp(*this); // to get around the fact that you can't resolve a const Data |
1945 |
temp.resolve(); |
1946 |
temp.saveDX(fileName); |
1947 |
return; |
1948 |
} |
1949 |
boost::python::dict args; |
1950 |
args["data"]=boost::python::object(this); |
1951 |
getDomain()->saveDX(fileName,args); |
1952 |
return; |
1953 |
} |
1954 |
|
1955 |
void |
1956 |
Data::saveVTK(std::string fileName) const |
1957 |
{ |
1958 |
if (isEmpty()) |
1959 |
{ |
1960 |
throw DataException("Error - Operations not permitted on instances of DataEmpty."); |
1961 |
} |
1962 |
if (isLazy()) |
1963 |
{ |
1964 |
Data temp(*this); // to get around the fact that you can't resolve a const Data |
1965 |
temp.resolve(); |
1966 |
temp.saveVTK(fileName); |
1967 |
return; |
1968 |
} |
1969 |
boost::python::dict args; |
1970 |
args["data"]=boost::python::object(this); |
1971 |
getDomain()->saveVTK(fileName,args); |
1972 |
return; |
1973 |
} |
1974 |
|
1975 |
Data& |
1976 |
Data::operator+=(const Data& right) |
1977 |
{ |
1978 |
if (isProtected()) { |
1979 |
throw DataException("Error - attempt to update protected Data object."); |
1980 |
} |
1981 |
if (isLazy() || right.isLazy()) |
1982 |
{ |
1983 |
DataLazy* c=new DataLazy(m_data,right.borrowDataPtr(),ADD); // for lazy + is equivalent to += |
1984 |
m_data=c->getPtr(); |
1985 |
return (*this); |
1986 |
} |
1987 |
else |
1988 |
{ |
1989 |
binaryOp(right,plus<double>()); |
1990 |
return (*this); |
1991 |
} |
1992 |
} |
1993 |
|
1994 |
Data& |
1995 |
Data::operator+=(const boost::python::object& right) |
1996 |
{ |
1997 |
if (isProtected()) { |
1998 |
throw DataException("Error - attempt to update protected Data object."); |
1999 |
} |
2000 |
Data tmp(right,getFunctionSpace(),false); |
2001 |
if (isLazy()) |
2002 |
{ |
2003 |
DataLazy* c=new DataLazy(m_data,tmp.borrowDataPtr(),ADD); // for lazy + is equivalent to += |
2004 |
m_data=c->getPtr(); |
2005 |
return (*this); |
2006 |
} |
2007 |
else |
2008 |
{ |
2009 |
binaryOp(tmp,plus<double>()); |
2010 |
return (*this); |
2011 |
} |
2012 |
} |
2013 |
|
2014 |
// Hmmm, operator= makes a deep copy but the copy constructor does not? |
2015 |
Data& |
2016 |
Data::operator=(const Data& other) |
2017 |
{ |
2018 |
copy(other); |
2019 |
return (*this); |
2020 |
} |
2021 |
|
2022 |
Data& |
2023 |
Data::operator-=(const Data& right) |
2024 |
{ |
2025 |
if (isProtected()) { |
2026 |
throw DataException("Error - attempt to update protected Data object."); |
2027 |
} |
2028 |
if (isLazy() || right.isLazy()) |
2029 |
{ |
2030 |
DataLazy* c=new DataLazy(m_data,right.borrowDataPtr(),SUB); // for lazy - is equivalent to -= |
2031 |
m_data=c->getPtr(); |
2032 |
return (*this); |
2033 |
} |
2034 |
else |
2035 |
{ |
2036 |
binaryOp(right,minus<double>()); |
2037 |
return (*this); |
2038 |
} |
2039 |
} |
2040 |
|
2041 |
Data& |
2042 |
Data::operator-=(const boost::python::object& right) |
2043 |
{ |
2044 |
if (isProtected()) { |
2045 |
throw DataException("Error - attempt to update protected Data object."); |
2046 |
} |
2047 |
Data tmp(right,getFunctionSpace(),false); |
2048 |
if (isLazy()) |
2049 |
{ |
2050 |
DataLazy* c=new DataLazy(m_data,tmp.borrowDataPtr(),SUB); // for lazy - is equivalent to -= |
2051 |
m_data=c->getPtr(); |
2052 |
return (*this); |
2053 |
} |
2054 |
else |
2055 |
{ |
2056 |
binaryOp(tmp,minus<double>()); |
2057 |
return (*this); |
2058 |
} |
2059 |
} |
2060 |
|
2061 |
Data& |
2062 |
Data::operator*=(const Data& right) |
2063 |
{ |
2064 |
if (isProtected()) { |
2065 |
throw DataException("Error - attempt to update protected Data object."); |
2066 |
} |
2067 |
if (isLazy() || right.isLazy()) |
2068 |
{ |
2069 |
DataLazy* c=new DataLazy(m_data,right.borrowDataPtr(),MUL); // for lazy * is equivalent to *= |
2070 |
m_data=c->getPtr(); |
2071 |
return (*this); |
2072 |
} |
2073 |
else |
2074 |
{ |
2075 |
binaryOp(right,multiplies<double>()); |
2076 |
return (*this); |
2077 |
} |
2078 |
} |
2079 |
|
2080 |
Data& |
2081 |
Data::operator*=(const boost::python::object& right) |
2082 |
{ |
2083 |
if (isProtected()) { |
2084 |
throw DataException("Error - attempt to update protected Data object."); |
2085 |
} |
2086 |
Data tmp(right,getFunctionSpace(),false); |
2087 |
if (isLazy()) |
2088 |
{ |
2089 |
DataLazy* c=new DataLazy(m_data,tmp.borrowDataPtr(),MUL); // for lazy * is equivalent to *= |
2090 |
m_data=c->getPtr(); |
2091 |
return (*this); |
2092 |
} |
2093 |
else |
2094 |
{ |
2095 |
binaryOp(tmp,multiplies<double>()); |
2096 |
return (*this); |
2097 |
} |
2098 |
} |
2099 |
|
2100 |
Data& |
2101 |
Data::operator/=(const Data& right) |
2102 |
{ |
2103 |
if (isProtected()) { |
2104 |
throw DataException("Error - attempt to update protected Data object."); |
2105 |
} |
2106 |
if (isLazy() || right.isLazy()) |
2107 |
{ |
2108 |
DataLazy* c=new DataLazy(m_data,right.borrowDataPtr(),DIV); // for lazy / is equivalent to /= |
2109 |
m_data=c->getPtr(); |
2110 |
return (*this); |
2111 |
} |
2112 |
else |
2113 |
{ |
2114 |
binaryOp(right,divides<double>()); |
2115 |
return (*this); |
2116 |
} |
2117 |
} |
2118 |
|
2119 |
Data& |
2120 |
Data::operator/=(const boost::python::object& right) |
2121 |
{ |
2122 |
if (isProtected()) { |
2123 |
throw DataException("Error - attempt to update protected Data object."); |
2124 |
} |
2125 |
Data tmp(right,getFunctionSpace(),false); |
2126 |
if (isLazy()) |
2127 |
{ |
2128 |
DataLazy* c=new DataLazy(m_data,tmp.borrowDataPtr(),DIV); // for lazy / is equivalent to /= |
2129 |
m_data=c->getPtr(); |
2130 |
return (*this); |
2131 |
} |
2132 |
else |
2133 |
{ |
2134 |
binaryOp(tmp,divides<double>()); |
2135 |
return (*this); |
2136 |
} |
2137 |
} |
2138 |
|
2139 |
Data |
2140 |
Data::rpowO(const boost::python::object& left) const |
2141 |
{ |
2142 |
Data left_d(left,*this); |
2143 |
return left_d.powD(*this); |
2144 |
} |
2145 |
|
2146 |
Data |
2147 |
Data::powO(const boost::python::object& right) const |
2148 |
{ |
2149 |
Data tmp(right,getFunctionSpace(),false); |
2150 |
return powD(tmp); |
2151 |
} |
2152 |
|
2153 |
Data |
2154 |
Data::powD(const Data& right) const |
2155 |
{ |
2156 |
if (isLazy() || right.isLazy()) |
2157 |
{ |
2158 |
DataLazy* c=new DataLazy(m_data,right.borrowDataPtr(),POW); |
2159 |
return Data(c); |
2160 |
} |
2161 |
return C_TensorBinaryOperation<double (*)(double, double)>(*this, right, ::pow); |
2162 |
} |
2163 |
|
2164 |
// |
2165 |
// NOTE: It is essential to specify the namespace this operator belongs to |
2166 |
Data |
2167 |
escript::operator+(const Data& left, const Data& right) |
2168 |
{ |
2169 |
if (left.isLazy() || right.isLazy()) |
2170 |
{ |
2171 |
DataLazy* c=new DataLazy(left.borrowDataPtr(),right.borrowDataPtr(),ADD); |
2172 |
return Data(c); |
2173 |
} |
2174 |
return C_TensorBinaryOperation(left, right, plus<double>()); |
2175 |
} |
2176 |
|
2177 |
// |
2178 |
// NOTE: It is essential to specify the namespace this operator belongs to |
2179 |
Data |
2180 |
escript::operator-(const Data& left, const Data& right) |
2181 |
{ |
2182 |
if (left.isLazy() || right.isLazy()) |
2183 |
{ |
2184 |
DataLazy* c=new DataLazy(left.borrowDataPtr(),right.borrowDataPtr(),SUB); |
2185 |
return Data(c); |
2186 |
} |
2187 |
return C_TensorBinaryOperation(left, right, minus<double>()); |
2188 |
} |
2189 |
|
2190 |
// |
2191 |
// NOTE: It is essential to specify the namespace this operator belongs to |
2192 |
Data |
2193 |
escript::operator*(const Data& left, const Data& right) |
2194 |
{ |
2195 |
if (left.isLazy() || right.isLazy()) |
2196 |
{ |
2197 |
DataLazy* c=new DataLazy(left.borrowDataPtr(),right.borrowDataPtr(),MUL); |
2198 |
return Data(c); |
2199 |
} |
2200 |
return C_TensorBinaryOperation(left, right, multiplies<double>()); |
2201 |
} |
2202 |
|
2203 |
// |
2204 |
// NOTE: It is essential to specify the namespace this operator belongs to |
2205 |
Data |
2206 |
escript::operator/(const Data& left, const Data& right) |
2207 |
{ |
2208 |
if (left.isLazy() || right.isLazy()) |
2209 |
{ |
2210 |
DataLazy* c=new DataLazy(left.borrowDataPtr(),right.borrowDataPtr(),DIV); |
2211 |
return Data(c); |
2212 |
} |
2213 |
return C_TensorBinaryOperation(left, right, divides<double>()); |
2214 |
} |
2215 |
|
2216 |
// |
2217 |
// NOTE: It is essential to specify the namespace this operator belongs to |
2218 |
Data |
2219 |
escript::operator+(const Data& left, const boost::python::object& right) |
2220 |
{ |
2221 |
if (left.isLazy()) |
2222 |
{ |
2223 |
DataLazy* c=new DataLazy(left.borrowDataPtr(),Data(right,left.getFunctionSpace(),false).borrowDataPtr(),ADD); |
2224 |
return Data(c); |
2225 |
} |
2226 |
return left+Data(right,left.getFunctionSpace(),false); |
2227 |
} |
2228 |
|
2229 |
// |
2230 |
// NOTE: It is essential to specify the namespace this operator belongs to |
2231 |
Data |
2232 |
escript::operator-(const Data& left, const boost::python::object& right) |
2233 |
{ |
2234 |
if (left.isLazy()) |
2235 |
{ |
2236 |
DataLazy* c=new DataLazy(left.borrowDataPtr(),Data(right,left.getFunctionSpace(),false).borrowDataPtr(),SUB); |
2237 |
return Data(c); |
2238 |
} |
2239 |
return left-Data(right,left.getFunctionSpace(),false); |
2240 |
} |
2241 |
|
2242 |
// |
2243 |
// NOTE: It is essential to specify the namespace this operator belongs to |
2244 |
Data |
2245 |
escript::operator*(const Data& left, const boost::python::object& right) |
2246 |
{ |
2247 |
if (left.isLazy()) |
2248 |
{ |
2249 |
DataLazy* c=new DataLazy(left.borrowDataPtr(),Data(right,left.getFunctionSpace(),false).borrowDataPtr(),MUL); |
2250 |
return Data(c); |
2251 |
} |
2252 |
return left*Data(right,left.getFunctionSpace(),false); |
2253 |
} |
2254 |
|
2255 |
// |
2256 |
// NOTE: It is essential to specify the namespace this operator belongs to |
2257 |
Data |
2258 |
escript::operator/(const Data& left, const boost::python::object& right) |
2259 |
{ |
2260 |
if (left.isLazy()) |
2261 |
{ |
2262 |
DataLazy* c=new DataLazy(left.borrowDataPtr(),Data(right,left.getFunctionSpace(),false).borrowDataPtr(),DIV); |
2263 |
return Data(c); |
2264 |
} |
2265 |
return left/Data(right,left.getFunctionSpace(),false); |
2266 |
} |
2267 |
|
2268 |
// |
2269 |
// NOTE: It is essential to specify the namespace this operator belongs to |
2270 |
Data |
2271 |
escript::operator+(const boost::python::object& left, const Data& right) |
2272 |
{ |
2273 |
if (right.isLazy()) |
2274 |
{ |
2275 |
DataLazy* c=new DataLazy(Data(left,right.getFunctionSpace(),false).borrowDataPtr(),right.borrowDataPtr(),ADD); |
2276 |
return Data(c); |
2277 |
} |
2278 |
return Data(left,right.getFunctionSpace(),false)+right; |
2279 |
} |
2280 |
|
2281 |
// |
2282 |
// NOTE: It is essential to specify the namespace this operator belongs to |
2283 |
Data |
2284 |
escript::operator-(const boost::python::object& left, const Data& right) |
2285 |
{ |
2286 |
if (right.isLazy()) |
2287 |
{ |
2288 |
DataLazy* c=new DataLazy(Data(left,right.getFunctionSpace(),false).borrowDataPtr(),right.borrowDataPtr(),SUB); |
2289 |
return Data(c); |
2290 |
} |
2291 |
return Data(left,right.getFunctionSpace(),false)-right; |
2292 |
} |
2293 |
|
2294 |
// |
2295 |
// NOTE: It is essential to specify the namespace this operator belongs to |
2296 |
Data |
2297 |
escript::operator*(const boost::python::object& left, const Data& right) |
2298 |
{ |
2299 |
if (right.isLazy()) |
2300 |
{ |
2301 |
DataLazy* c=new DataLazy(Data(left,right.getFunctionSpace(),false).borrowDataPtr(),right.borrowDataPtr(),MUL); |
2302 |
return Data(c); |
2303 |
} |
2304 |
return Data(left,right.getFunctionSpace(),false)*right; |
2305 |
} |
2306 |
|
2307 |
// |
2308 |
// NOTE: It is essential to specify the namespace this operator belongs to |
2309 |
Data |
2310 |
escript::operator/(const boost::python::object& left, const Data& right) |
2311 |
{ |
2312 |
if (right.isLazy()) |
2313 |
{ |
2314 |
DataLazy* c=new DataLazy(Data(left,right.getFunctionSpace(),false).borrowDataPtr(),right.borrowDataPtr(),DIV); |
2315 |
return Data(c); |
2316 |
} |
2317 |
return Data(left,right.getFunctionSpace(),false)/right; |
2318 |
} |
2319 |
|
2320 |
|
2321 |
/* TODO */ |
2322 |
/* global reduction */ |
2323 |
Data |
2324 |
Data::getItem(const boost::python::object& key) const |
2325 |
{ |
2326 |
|
2327 |
DataTypes::RegionType slice_region=DataTypes::getSliceRegion(getDataPointShape(),key); |
2328 |
|
2329 |
if (slice_region.size()!=getDataPointRank()) { |
2330 |
throw DataException("Error - slice size does not match Data rank."); |
2331 |
} |
2332 |
|
2333 |
return getSlice(slice_region); |
2334 |
} |
2335 |
|
2336 |
/* TODO */ |
2337 |
/* global reduction */ |
2338 |
Data |
2339 |
Data::getSlice(const DataTypes::RegionType& region) const |
2340 |
{ |
2341 |
return Data(*this,region); |
2342 |
} |
2343 |
|
2344 |
/* TODO */ |
2345 |
/* global reduction */ |
2346 |
void |
2347 |
Data::setItemO(const boost::python::object& key, |
2348 |
const boost::python::object& value) |
2349 |
{ |
2350 |
Data tempData(value,getFunctionSpace()); |
2351 |
setItemD(key,tempData); |
2352 |
} |
2353 |
|
2354 |
void |
2355 |
Data::setItemD(const boost::python::object& key, |
2356 |
const Data& value) |
2357 |
{ |
2358 |
// const DataArrayView& view=getPointDataView(); |
2359 |
|
2360 |
DataTypes::RegionType slice_region=DataTypes::getSliceRegion(getDataPointShape(),key); |
2361 |
if (slice_region.size()!=getDataPointRank()) { |
2362 |
throw DataException("Error - slice size does not match Data rank."); |
2363 |
} |
2364 |
if (getFunctionSpace()!=value.getFunctionSpace()) { |
2365 |
setSlice(Data(value,getFunctionSpace()),slice_region); |
2366 |
} else { |
2367 |
setSlice(value,slice_region); |
2368 |
} |
2369 |
} |
2370 |
|
2371 |
void |
2372 |
Data::setSlice(const Data& value, |
2373 |
const DataTypes::RegionType& region) |
2374 |
{ |
2375 |
if (isProtected()) { |
2376 |
throw DataException("Error - attempt to update protected Data object."); |
2377 |
} |
2378 |
FORCERESOLVE; |
2379 |
/* if (isLazy()) |
2380 |
{ |
2381 |
throw DataException("Error - setSlice not permitted on lazy data."); |
2382 |
}*/ |
2383 |
Data tempValue(value); |
2384 |
typeMatchLeft(tempValue); |
2385 |
typeMatchRight(tempValue); |
2386 |
getReady()->setSlice(tempValue.m_data.get(),region); |
2387 |
} |
2388 |
|
2389 |
void |
2390 |
Data::typeMatchLeft(Data& right) const |
2391 |
{ |
2392 |
if (right.isLazy() && !isLazy()) |
2393 |
{ |
2394 |
right.resolve(); |
2395 |
} |
2396 |
if (isExpanded()){ |
2397 |
right.expand(); |
2398 |
} else if (isTagged()) { |
2399 |
if (right.isConstant()) { |
2400 |
right.tag(); |
2401 |
} |
2402 |
} |
2403 |
} |
2404 |
|
2405 |
void |
2406 |
Data::typeMatchRight(const Data& right) |
2407 |
{ |
2408 |
if (isLazy() && !right.isLazy()) |
2409 |
{ |
2410 |
resolve(); |
2411 |
} |
2412 |
if (isTagged()) { |
2413 |
if (right.isExpanded()) { |
2414 |
expand(); |
2415 |
} |
2416 |
} else if (isConstant()) { |
2417 |
if (right.isExpanded()) { |
2418 |
expand(); |
2419 |
} else if (right.isTagged()) { |
2420 |
tag(); |
2421 |
} |
2422 |
} |
2423 |
} |
2424 |
|
2425 |
void |
2426 |
Data::setTaggedValueByName(std::string name, |
2427 |
const boost::python::object& value) |
2428 |
{ |
2429 |
if (getFunctionSpace().getDomain()->isValidTagName(name)) { |
2430 |
FORCERESOLVE; |
2431 |
int tagKey=getFunctionSpace().getDomain()->getTag(name); |
2432 |
setTaggedValue(tagKey,value); |
2433 |
} |
2434 |
} |
2435 |
void |
2436 |
Data::setTaggedValue(int tagKey, |
2437 |
const boost::python::object& value) |
2438 |
{ |
2439 |
if (isProtected()) { |
2440 |
throw DataException("Error - attempt to update protected Data object."); |
2441 |
} |
2442 |
// |
2443 |
// Ensure underlying data object is of type DataTagged |
2444 |
FORCERESOLVE; |
2445 |
if (isConstant()) tag(); |
2446 |
numeric::array asNumArray(value); |
2447 |
|
2448 |
// extract the shape of the numarray |
2449 |
DataTypes::ShapeType tempShape; |
2450 |
for (int i=0; i < asNumArray.getrank(); i++) { |
2451 |
tempShape.push_back(extract<int>(asNumArray.getshape()[i])); |
2452 |
} |
2453 |
|
2454 |
DataVector temp_data2; |
2455 |
temp_data2.copyFromNumArray(asNumArray); |
2456 |
|
2457 |
m_data->setTaggedValue(tagKey,tempShape, temp_data2); |
2458 |
} |
2459 |
|
2460 |
|
2461 |
void |
2462 |
Data::setTaggedValueFromCPP(int tagKey, |
2463 |
const DataTypes::ShapeType& pointshape, |
2464 |
const DataTypes::ValueType& value, |
2465 |
int dataOffset) |
2466 |
{ |
2467 |
if (isProtected()) { |
2468 |
throw DataException("Error - attempt to update protected Data object."); |
2469 |
} |
2470 |
// |
2471 |
// Ensure underlying data object is of type DataTagged |
2472 |
FORCERESOLVE; |
2473 |
if (isConstant()) tag(); |
2474 |
// |
2475 |
// Call DataAbstract::setTaggedValue |
2476 |
m_data->setTaggedValue(tagKey,pointshape, value, dataOffset); |
2477 |
} |
2478 |
|
2479 |
int |
2480 |
Data::getTagNumber(int dpno) |
2481 |
{ |
2482 |
if (isEmpty()) |
2483 |
{ |
2484 |
throw DataException("Error - operation not permitted on instances of DataEmpty."); |
2485 |
} |
2486 |
return getFunctionSpace().getTagFromDataPointNo(dpno); |
2487 |
} |
2488 |
|
2489 |
|
2490 |
ostream& escript::operator<<(ostream& o, const Data& data) |
2491 |
{ |
2492 |
o << data.toString(); |
2493 |
return o; |
2494 |
} |
2495 |
|
2496 |
Data |
2497 |
escript::C_GeneralTensorProduct(Data& arg_0, |
2498 |
Data& arg_1, |
2499 |
int axis_offset, |
2500 |
int transpose) |
2501 |
{ |
2502 |
// General tensor product: res(SL x SR) = arg_0(SL x SM) * arg_1(SM x SR) |
2503 |
// SM is the product of the last axis_offset entries in arg_0.getShape(). |
2504 |
|
2505 |
// deal with any lazy data |
2506 |
// if (arg_0.isLazy()) {arg_0.resolve();} |
2507 |
// if (arg_1.isLazy()) {arg_1.resolve();} |
2508 |
if (arg_0.isLazy() || arg_1.isLazy()) |
2509 |
{ |
2510 |
DataLazy* c=new DataLazy(arg_0.borrowDataPtr(), arg_1.borrowDataPtr(), PROD, axis_offset,transpose); |
2511 |
return Data(c); |
2512 |
} |
2513 |
|
2514 |
// Interpolate if necessary and find an appropriate function space |
2515 |
Data arg_0_Z, arg_1_Z; |
2516 |
if (arg_0.getFunctionSpace()!=arg_1.getFunctionSpace()) { |
2517 |
if (arg_0.probeInterpolation(arg_1.getFunctionSpace())) { |
2518 |
arg_0_Z = arg_0.interpolate(arg_1.getFunctionSpace()); |
2519 |
arg_1_Z = Data(arg_1); |
2520 |
} |
2521 |
else if (arg_1.probeInterpolation(arg_0.getFunctionSpace())) { |
2522 |
arg_1_Z=arg_1.interpolate(arg_0.getFunctionSpace()); |
2523 |
arg_0_Z =Data(arg_0); |
2524 |
} |
2525 |
else { |
2526 |
throw DataException("Error - C_GeneralTensorProduct: arguments have incompatible function spaces."); |
2527 |
} |
2528 |
} else { |
2529 |
arg_0_Z = Data(arg_0); |
2530 |
arg_1_Z = Data(arg_1); |
2531 |
} |
2532 |
// Get rank and shape of inputs |
2533 |
int rank0 = arg_0_Z.getDataPointRank(); |
2534 |
int rank1 = arg_1_Z.getDataPointRank(); |
2535 |
const DataTypes::ShapeType& shape0 = arg_0_Z.getDataPointShape(); |
2536 |
const DataTypes::ShapeType& shape1 = arg_1_Z.getDataPointShape(); |
2537 |
|
2538 |
// Prepare for the loops of the product and verify compatibility of shapes |
2539 |
int start0=0, start1=0; |
2540 |
if (transpose == 0) {} |
2541 |
else if (transpose == 1) { start0 = axis_offset; } |
2542 |
else if (transpose == 2) { start1 = rank1-axis_offset; } |
2543 |
else { throw DataException("C_GeneralTensorProduct: Error - transpose should be 0, 1 or 2"); } |
2544 |
|
2545 |
|
2546 |
// Adjust the shapes for transpose |
2547 |
DataTypes::ShapeType tmpShape0(rank0); // pre-sizing the vectors rather |
2548 |
DataTypes::ShapeType tmpShape1(rank1); // than using push_back |
2549 |
for (int i=0; i<rank0; i++) { tmpShape0[i]=shape0[(i+start0)%rank0]; } |
2550 |
for (int i=0; i<rank1; i++) { tmpShape1[i]=shape1[(i+start1)%rank1]; } |
2551 |
|
2552 |
#if 0 |
2553 |
// For debugging: show shape after transpose |
2554 |
char tmp[100]; |
2555 |
std::string shapeStr; |
2556 |
shapeStr = "("; |
2557 |
for (int i=0; i<rank0; i++) { sprintf(tmp, "%d,", tmpShape0[i]); shapeStr += tmp; } |
2558 |
shapeStr += ")"; |
2559 |
cout << "C_GeneralTensorProduct: Shape of arg0 is " << shapeStr << endl; |
2560 |
shapeStr = "("; |
2561 |
for (int i=0; i<rank1; i++) { sprintf(tmp, "%d,", tmpShape1[i]); shapeStr += tmp; } |
2562 |
shapeStr += ")"; |
2563 |
cout << "C_GeneralTensorProduct: Shape of arg1 is " << shapeStr << endl; |
2564 |
#endif |
2565 |
|
2566 |
// Prepare for the loops of the product |
2567 |
int SL=1, SM=1, SR=1; |
2568 |
for (int i=0; i<rank0-axis_offset; i++) { |
2569 |
SL *= tmpShape0[i]; |
2570 |
} |
2571 |
for (int i=rank0-axis_offset; i<rank0; i++) { |
2572 |
if (tmpShape0[i] != tmpShape1[i-(rank0-axis_offset)]) { |
2573 |
throw DataException("C_GeneralTensorProduct: Error - incompatible shapes"); |
2574 |
} |
2575 |
SM *= tmpShape0[i]; |
2576 |
} |
2577 |
for (int i=axis_offset; i<rank1; i++) { |
2578 |
SR *= tmpShape1[i]; |
2579 |
} |
2580 |
|
2581 |
// Define the shape of the output (rank of shape is the sum of the loop ranges below) |
2582 |
DataTypes::ShapeType shape2(rank0+rank1-2*axis_offset); |
2583 |
{ // block to limit the scope of out_index |
2584 |
int out_index=0; |
2585 |
for (int i=0; i<rank0-axis_offset; i++, ++out_index) { shape2[out_index]=tmpShape0[i]; } // First part of arg_0_Z |
2586 |
for (int i=axis_offset; i<rank1; i++, ++out_index) { shape2[out_index]=tmpShape1[i]; } // Last part of arg_1_Z |
2587 |
} |
2588 |
|
2589 |
if (shape2.size()>ESCRIPT_MAX_DATA_RANK) |
2590 |
{ |
2591 |
ostringstream os; |
2592 |
os << "C_GeneralTensorProduct: Error - Attempt to create a rank " << shape2.size() << " object. The maximum rank is " << ESCRIPT_MAX_DATA_RANK << "."; |
2593 |
throw DataException(os.str()); |
2594 |
} |
2595 |
|
2596 |
// Declare output Data object |
2597 |
Data res; |
2598 |
|
2599 |
if (arg_0_Z.isConstant() && arg_1_Z.isConstant()) { |
2600 |
res = Data(0.0, shape2, arg_1_Z.getFunctionSpace()); // DataConstant output |
2601 |
double *ptr_0 = &(arg_0_Z.getDataAtOffset(0)); |
2602 |
double *ptr_1 = &(arg_1_Z.getDataAtOffset(0)); |
2603 |
double *ptr_2 = &(res.getDataAtOffset(0)); |
2604 |
matrix_matrix_product(SL, SM, SR, ptr_0, ptr_1, ptr_2, transpose); |
2605 |
} |
2606 |
else if (arg_0_Z.isConstant() && arg_1_Z.isTagged()) { |
2607 |
|
2608 |
// Prepare the DataConstant input |
2609 |
DataConstant* tmp_0=dynamic_cast<DataConstant*>(arg_0_Z.borrowData()); |
2610 |
if (tmp_0==0) { throw DataException("GTP Programming error - casting to DataConstant."); } |
2611 |
|
2612 |
// Borrow DataTagged input from Data object |
2613 |
DataTagged* tmp_1=dynamic_cast<DataTagged*>(arg_1_Z.borrowData()); |
2614 |
if (tmp_1==0) { throw DataException("GTP_1 Programming error - casting to DataTagged."); } |
2615 |
|
2616 |
// Prepare a DataTagged output 2 |
2617 |
res = Data(0.0, shape2, arg_1_Z.getFunctionSpace()); // DataTagged output |
2618 |
res.tag(); |
2619 |
DataTagged* tmp_2=dynamic_cast<DataTagged*>(res.borrowData()); |
2620 |
if (tmp_2==0) { throw DataException("GTP Programming error - casting to DataTagged."); } |
2621 |
|
2622 |
// Prepare offset into DataConstant |
2623 |
int offset_0 = tmp_0->getPointOffset(0,0); |
2624 |
double *ptr_0 = &(arg_0_Z.getDataAtOffset(offset_0)); |
2625 |
// Get the views |
2626 |
// DataArrayView view_1 = tmp_1->getDefaultValue(); |
2627 |
// DataArrayView view_2 = tmp_2->getDefaultValue(); |
2628 |
// // Get the pointers to the actual data |
2629 |
// double *ptr_1 = &((view_1.getData())[0]); |
2630 |
// double *ptr_2 = &((view_2.getData())[0]); |
2631 |
|
2632 |
double *ptr_1 = &(tmp_1->getDefaultValue(0)); |
2633 |
double *ptr_2 = &(tmp_2->getDefaultValue(0)); |
2634 |
|
2635 |
|
2636 |
// Compute an MVP for the default |
2637 |
matrix_matrix_product(SL, SM, SR, ptr_0, ptr_1, ptr_2, transpose); |
2638 |
// Compute an MVP for each tag |
2639 |
const DataTagged::DataMapType& lookup_1=tmp_1->getTagLookup(); |
2640 |
DataTagged::DataMapType::const_iterator i; // i->first is a tag, i->second is an offset into memory |
2641 |
for (i=lookup_1.begin();i!=lookup_1.end();i++) { |
2642 |
tmp_2->addTag(i->first); |
2643 |
// DataArrayView view_1 = tmp_1->getDataPointByTag(i->first); |
2644 |
// DataArrayView view_2 = tmp_2->getDataPointByTag(i->first); |
2645 |
// double *ptr_1 = &view_1.getData(0); |
2646 |
// double *ptr_2 = &view_2.getData(0); |
2647 |
|
2648 |
double *ptr_1 = &(tmp_1->getDataByTag(i->first,0)); |
2649 |
double *ptr_2 = &(tmp_2->getDataByTag(i->first,0)); |
2650 |
|
2651 |
matrix_matrix_product(SL, SM, SR, ptr_0, ptr_1, ptr_2, transpose); |
2652 |
} |
2653 |
|
2654 |
} |
2655 |
else if (arg_0_Z.isConstant() && arg_1_Z.isExpanded()) { |
2656 |
|
2657 |
res = Data(0.0, shape2, arg_1_Z.getFunctionSpace(),true); // DataExpanded output |
2658 |
DataConstant* tmp_0=dynamic_cast<DataConstant*>(arg_0_Z.borrowData()); |
2659 |
DataExpanded* tmp_1=dynamic_cast<DataExpanded*>(arg_1_Z.borrowData()); |
2660 |
DataExpanded* tmp_2=dynamic_cast<DataExpanded*>(res.borrowData()); |
2661 |
if (tmp_0==0) { throw DataException("GTP Programming error - casting to DataConstant."); } |
2662 |
if (tmp_1==0) { throw DataException("GTP Programming error - casting to DataExpanded."); } |
2663 |
if (tmp_2==0) { throw DataException("GTP Programming error - casting to DataExpanded."); } |
2664 |
int sampleNo_1,dataPointNo_1; |
2665 |
int numSamples_1 = arg_1_Z.getNumSamples(); |
2666 |
int numDataPointsPerSample_1 = arg_1_Z.getNumDataPointsPerSample(); |
2667 |
int offset_0 = tmp_0->getPointOffset(0,0); |
2668 |
#pragma omp parallel for private(sampleNo_1,dataPointNo_1) schedule(static) |
2669 |
for (sampleNo_1 = 0; sampleNo_1 < numSamples_1; sampleNo_1++) { |
2670 |
for (dataPointNo_1 = 0; dataPointNo_1 < numDataPointsPerSample_1; dataPointNo_1++) { |
2671 |
int offset_1 = tmp_1->getPointOffset(sampleNo_1,dataPointNo_1); |
2672 |
int offset_2 = tmp_2->getPointOffset(sampleNo_1,dataPointNo_1); |
2673 |
double *ptr_0 = &(arg_0_Z.getDataAtOffset(offset_0)); |
2674 |
double *ptr_1 = &(arg_1_Z.getDataAtOffset(offset_1)); |
2675 |
double *ptr_2 = &(res.getDataAtOffset(offset_2)); |
2676 |
matrix_matrix_product(SL, SM, SR, ptr_0, ptr_1, ptr_2, transpose); |
2677 |
} |
2678 |
} |
2679 |
|
2680 |
} |
2681 |
else if (arg_0_Z.isTagged() && arg_1_Z.isConstant()) { |
2682 |
|
2683 |
// Borrow DataTagged input from Data object |
2684 |
DataTagged* tmp_0=dynamic_cast<DataTagged*>(arg_0_Z.borrowData()); |
2685 |
if (tmp_0==0) { throw DataException("GTP_0 Programming error - casting to DataTagged."); } |
2686 |
|
2687 |
// Prepare the DataConstant input |
2688 |
DataConstant* tmp_1=dynamic_cast<DataConstant*>(arg_1_Z.borrowData()); |
2689 |
if (tmp_1==0) { throw DataException("GTP Programming error - casting to DataConstant."); } |
2690 |
|
2691 |
// Prepare a DataTagged output 2 |
2692 |
res = Data(0.0, shape2, arg_0_Z.getFunctionSpace()); // DataTagged output |
2693 |
res.tag(); |
2694 |
DataTagged* tmp_2=dynamic_cast<DataTagged*>(res.borrowData()); |
2695 |
if (tmp_2==0) { throw DataException("GTP Programming error - casting to DataTagged."); } |
2696 |
|
2697 |
// Prepare offset into DataConstant |
2698 |
int offset_1 = tmp_1->getPointOffset(0,0); |
2699 |
double *ptr_1 = &(arg_1_Z.getDataAtOffset(offset_1)); |
2700 |
// Get the views |
2701 |
// DataArrayView view_0 = tmp_0->getDefaultValue(); |
2702 |
// DataArrayView view_2 = tmp_2->getDefaultValue(); |
2703 |
// // Get the pointers to the actual data |
2704 |
// double *ptr_0 = &((view_0.getData())[0]); |
2705 |
// double *ptr_2 = &((view_2.getData())[0]); |
2706 |
|
2707 |
double *ptr_0 = &(tmp_0->getDefaultValue(0)); |
2708 |
double *ptr_2 = &(tmp_2->getDefaultValue(0)); |
2709 |
|
2710 |
// Compute an MVP for the default |
2711 |
matrix_matrix_product(SL, SM, SR, ptr_0, ptr_1, ptr_2, transpose); |
2712 |
// Compute an MVP for each tag |
2713 |
const DataTagged::DataMapType& lookup_0=tmp_0->getTagLookup(); |
2714 |
DataTagged::DataMapType::const_iterator i; // i->first is a tag, i->second is an offset into memory |
2715 |
for (i=lookup_0.begin();i!=lookup_0.end();i++) { |
2716 |
// tmp_2->addTaggedValue(i->first,tmp_2->getDefaultValue()); |
2717 |
// DataArrayView view_0 = tmp_0->getDataPointByTag(i->first); |
2718 |
// DataArrayView view_2 = tmp_2->getDataPointByTag(i->first); |
2719 |
// double *ptr_0 = &view_0.getData(0); |
2720 |
// double *ptr_2 = &view_2.getData(0); |
2721 |
|
2722 |
tmp_2->addTag(i->first); |
2723 |
double *ptr_0 = &(tmp_0->getDataByTag(i->first,0)); |
2724 |
double *ptr_2 = &(tmp_2->getDataByTag(i->first,0)); |
2725 |
matrix_matrix_product(SL, SM, SR, ptr_0, ptr_1, ptr_2, transpose); |
2726 |
} |
2727 |
|
2728 |
} |
2729 |
else if (arg_0_Z.isTagged() && arg_1_Z.isTagged()) { |
2730 |
|
2731 |
// Borrow DataTagged input from Data object |
2732 |
DataTagged* tmp_0=dynamic_cast<DataTagged*>(arg_0_Z.borrowData()); |
2733 |
if (tmp_0==0) { throw DataException("GTP Programming error - casting to DataTagged."); } |
2734 |
|
2735 |
// Borrow DataTagged input from Data object |
2736 |
DataTagged* tmp_1=dynamic_cast<DataTagged*>(arg_1_Z.borrowData()); |
2737 |
if (tmp_1==0) { throw DataException("GTP Programming error - casting to DataTagged."); } |
2738 |
|
2739 |
// Prepare a DataTagged output 2 |
2740 |
res = Data(0.0, shape2, arg_1_Z.getFunctionSpace()); |
2741 |
res.tag(); // DataTagged output |
2742 |
DataTagged* tmp_2=dynamic_cast<DataTagged*>(res.borrowData()); |
2743 |
if (tmp_2==0) { throw DataException("GTP Programming error - casting to DataTagged."); } |
2744 |
|
2745 |
// // Get the views |
2746 |
// DataArrayView view_0 = tmp_0->getDefaultValue(); |
2747 |
// DataArrayView view_1 = tmp_1->getDefaultValue(); |
2748 |
// DataArrayView view_2 = tmp_2->getDefaultValue(); |
2749 |
// // Get the pointers to the actual data |
2750 |
// double *ptr_0 = &((view_0.getData())[0]); |
2751 |
// double *ptr_1 = &((view_1.getData())[0]); |
2752 |
// double *ptr_2 = &((view_2.getData())[0]); |
2753 |
|
2754 |
double *ptr_0 = &(tmp_0->getDefaultValue(0)); |
2755 |
double *ptr_1 = &(tmp_1->getDefaultValue(0)); |
2756 |
double *ptr_2 = &(tmp_2->getDefaultValue(0)); |
2757 |
|
2758 |
|
2759 |
// Compute an MVP for the default |
2760 |
matrix_matrix_product(SL, SM, SR, ptr_0, ptr_1, ptr_2, transpose); |
2761 |
// Merge the tags |
2762 |
DataTagged::DataMapType::const_iterator i; // i->first is a tag, i->second is an offset into memory |
2763 |
const DataTagged::DataMapType& lookup_0=tmp_0->getTagLookup(); |
2764 |
const DataTagged::DataMapType& lookup_1=tmp_1->getTagLookup(); |
2765 |
for (i=lookup_0.begin();i!=lookup_0.end();i++) { |
2766 |
tmp_2->addTag(i->first); // use tmp_2 to get correct shape |
2767 |
} |
2768 |
for (i=lookup_1.begin();i!=lookup_1.end();i++) { |
2769 |
tmp_2->addTag(i->first); |
2770 |
} |
2771 |
// Compute an MVP for each tag |
2772 |
const DataTagged::DataMapType& lookup_2=tmp_2->getTagLookup(); |
2773 |
for (i=lookup_2.begin();i!=lookup_2.end();i++) { |
2774 |
// DataArrayView view_0 = tmp_0->getDataPointByTag(i->first); |
2775 |
// DataArrayView view_1 = tmp_1->getDataPointByTag(i->first); |
2776 |
// DataArrayView view_2 = tmp_2->getDataPointByTag(i->first); |
2777 |
// double *ptr_0 = &view_0.getData(0); |
2778 |
// double *ptr_1 = &view_1.getData(0); |
2779 |
// double *ptr_2 = &view_2.getData(0); |
2780 |
|
2781 |
double *ptr_0 = &(tmp_0->getDataByTag(i->first,0)); |
2782 |
double *ptr_1 = &(tmp_1->getDataByTag(i->first,0)); |
2783 |
double *ptr_2 = &(tmp_2->getDataByTag(i->first,0)); |
2784 |
|
2785 |
matrix_matrix_product(SL, SM, SR, ptr_0, ptr_1, ptr_2, transpose); |
2786 |
} |
2787 |
|
2788 |
} |
2789 |
else if (arg_0_Z.isTagged() && arg_1_Z.isExpanded()) { |
2790 |
|
2791 |
// After finding a common function space above the two inputs have the same numSamples and num DPPS |
2792 |
res = Data(0.0, shape2, arg_1_Z.getFunctionSpace(),true); // DataExpanded output |
2793 |
DataTagged* tmp_0=dynamic_cast<DataTagged*>(arg_0_Z.borrowData()); |
2794 |
DataExpanded* tmp_1=dynamic_cast<DataExpanded*>(arg_1_Z.borrowData()); |
2795 |
DataExpanded* tmp_2=dynamic_cast<DataExpanded*>(res.borrowData()); |
2796 |
if (tmp_0==0) { throw DataException("GTP Programming error - casting to DataTagged."); } |
2797 |
if (tmp_1==0) { throw DataException("GTP Programming error - casting to DataExpanded."); } |
2798 |
if (tmp_2==0) { throw DataException("GTP Programming error - casting to DataExpanded."); } |
2799 |
int sampleNo_0,dataPointNo_0; |
2800 |
int numSamples_0 = arg_0_Z.getNumSamples(); |
2801 |
int numDataPointsPerSample_0 = arg_0_Z.getNumDataPointsPerSample(); |
2802 |
#pragma omp parallel for private(sampleNo_0,dataPointNo_0) schedule(static) |
2803 |
for (sampleNo_0 = 0; sampleNo_0 < numSamples_0; sampleNo_0++) { |
2804 |
int offset_0 = tmp_0->getPointOffset(sampleNo_0,0); // They're all the same, so just use #0 |
2805 |
double *ptr_0 = &(arg_0_Z.getDataAtOffset(offset_0)); |
2806 |
for (dataPointNo_0 = 0; dataPointNo_0 < numDataPointsPerSample_0; dataPointNo_0++) { |
2807 |
int offset_1 = tmp_1->getPointOffset(sampleNo_0,dataPointNo_0); |
2808 |
int offset_2 = tmp_2->getPointOffset(sampleNo_0,dataPointNo_0); |
2809 |
double *ptr_1 = &(arg_1_Z.getDataAtOffset(offset_1)); |
2810 |
double *ptr_2 = &(res.getDataAtOffset(offset_2)); |
2811 |
matrix_matrix_product(SL, SM, SR, ptr_0, ptr_1, ptr_2, transpose); |
2812 |
} |
2813 |
} |
2814 |
|
2815 |
} |
2816 |
else if (arg_0_Z.isExpanded() && arg_1_Z.isConstant()) { |
2817 |
|
2818 |
res = Data(0.0, shape2, arg_1_Z.getFunctionSpace(),true); // DataExpanded output |
2819 |
DataExpanded* tmp_0=dynamic_cast<DataExpanded*>(arg_0_Z.borrowData()); |
2820 |
DataConstant* tmp_1=dynamic_cast<DataConstant*>(arg_1_Z.borrowData()); |
2821 |
DataExpanded* tmp_2=dynamic_cast<DataExpanded*>(res.borrowData()); |
2822 |
if (tmp_0==0) { throw DataException("GTP Programming error - casting to DataExpanded."); } |
2823 |
if (tmp_1==0) { throw DataException("GTP Programming error - casting to DataConstant."); } |
2824 |
if (tmp_2==0) { throw DataException("GTP Programming error - casting to DataExpanded."); } |
2825 |
int sampleNo_0,dataPointNo_0; |
2826 |
int numSamples_0 = arg_0_Z.getNumSamples(); |
2827 |
int numDataPointsPerSample_0 = arg_0_Z.getNumDataPointsPerSample(); |
2828 |
int offset_1 = tmp_1->getPointOffset(0,0); |
2829 |
#pragma omp parallel for private(sampleNo_0,dataPointNo_0) schedule(static) |
2830 |
for (sampleNo_0 = 0; sampleNo_0 < numSamples_0; sampleNo_0++) { |
2831 |
for (dataPointNo_0 = 0; dataPointNo_0 < numDataPointsPerSample_0; dataPointNo_0++) { |
2832 |
int offset_0 = tmp_0->getPointOffset(sampleNo_0,dataPointNo_0); |
2833 |
int offset_2 = tmp_2->getPointOffset(sampleNo_0,dataPointNo_0); |
2834 |
double *ptr_0 = &(arg_0_Z.getDataAtOffset(offset_0)); |
2835 |
double *ptr_1 = &(arg_1_Z.getDataAtOffset(offset_1)); |
2836 |
double *ptr_2 = &(res.getDataAtOffset(offset_2)); |
2837 |
matrix_matrix_product(SL, SM, SR, ptr_0, ptr_1, ptr_2, transpose); |
2838 |
} |
2839 |
} |
2840 |
|
2841 |
|
2842 |
} |
2843 |
else if (arg_0_Z.isExpanded() && arg_1_Z.isTagged()) { |
2844 |
|
2845 |
// After finding a common function space above the two inputs have the same numSamples and num DPPS |
2846 |
res = Data(0.0, shape2, arg_1_Z.getFunctionSpace(),true); // DataExpanded output |
2847 |
DataExpanded* tmp_0=dynamic_cast<DataExpanded*>(arg_0_Z.borrowData()); |
2848 |
DataTagged* tmp_1=dynamic_cast<DataTagged*>(arg_1_Z.borrowData()); |
2849 |
DataExpanded* tmp_2=dynamic_cast<DataExpanded*>(res.borrowData()); |
2850 |
if (tmp_0==0) { throw DataException("GTP Programming error - casting to DataExpanded."); } |
2851 |
if (tmp_1==0) { throw DataException("GTP Programming error - casting to DataTagged."); } |
2852 |
if (tmp_2==0) { throw DataException("GTP Programming error - casting to DataExpanded."); } |
2853 |
int sampleNo_0,dataPointNo_0; |
2854 |
int numSamples_0 = arg_0_Z.getNumSamples(); |
2855 |
int numDataPointsPerSample_0 = arg_0_Z.getNumDataPointsPerSample(); |
2856 |
#pragma omp parallel for private(sampleNo_0,dataPointNo_0) schedule(static) |
2857 |
for (sampleNo_0 = 0; sampleNo_0 < numSamples_0; sampleNo_0++) { |
2858 |
int offset_1 = tmp_1->getPointOffset(sampleNo_0,0); |
2859 |
double *ptr_1 = &(arg_1_Z.getDataAtOffset(offset_1)); |
2860 |
for (dataPointNo_0 = 0; dataPointNo_0 < numDataPointsPerSample_0; dataPointNo_0++) { |
2861 |
int offset_0 = tmp_0->getPointOffset(sampleNo_0,dataPointNo_0); |
2862 |
int offset_2 = tmp_2->getPointOffset(sampleNo_0,dataPointNo_0); |
2863 |
double *ptr_0 = &(arg_0_Z.getDataAtOffset(offset_0)); |
2864 |
double *ptr_2 = &(res.getDataAtOffset(offset_2)); |
2865 |
matrix_matrix_product(SL, SM, SR, ptr_0, ptr_1, ptr_2, transpose); |
2866 |
} |
2867 |
} |
2868 |
|
2869 |
} |
2870 |
else if (arg_0_Z.isExpanded() && arg_1_Z.isExpanded()) { |
2871 |
|
2872 |
// After finding a common function space above the two inputs have the same numSamples and num DPPS |
2873 |
res = Data(0.0, shape2, arg_1_Z.getFunctionSpace(),true); // DataExpanded output |
2874 |
DataExpanded* tmp_0=dynamic_cast<DataExpanded*>(arg_0_Z.borrowData()); |
2875 |
DataExpanded* tmp_1=dynamic_cast<DataExpanded*>(arg_1_Z.borrowData()); |
2876 |
DataExpanded* tmp_2=dynamic_cast<DataExpanded*>(res.borrowData()); |
2877 |
if (tmp_0==0) { throw DataException("GTP Programming error - casting to DataExpanded."); } |
2878 |
if (tmp_1==0) { throw DataException("GTP Programming error - casting to DataExpanded."); } |
2879 |
if (tmp_2==0) { throw DataException("GTP Programming error - casting to DataExpanded."); } |
2880 |
int sampleNo_0,dataPointNo_0; |
2881 |
int numSamples_0 = arg_0_Z.getNumSamples(); |
2882 |
int numDataPointsPerSample_0 = arg_0_Z.getNumDataPointsPerSample(); |
2883 |
#pragma omp parallel for private(sampleNo_0,dataPointNo_0) schedule(static) |
2884 |
for (sampleNo_0 = 0; sampleNo_0 < numSamples_0; sampleNo_0++) { |
2885 |
for (dataPointNo_0 = 0; dataPointNo_0 < numDataPointsPerSample_0; dataPointNo_0++) { |
2886 |
int offset_0 = tmp_0->getPointOffset(sampleNo_0,dataPointNo_0); |
2887 |
int offset_1 = tmp_1->getPointOffset(sampleNo_0,dataPointNo_0); |
2888 |
int offset_2 = tmp_2->getPointOffset(sampleNo_0,dataPointNo_0); |
2889 |
double *ptr_0 = &(arg_0_Z.getDataAtOffset(offset_0)); |
2890 |
double *ptr_1 = &(arg_1_Z.getDataAtOffset(offset_1)); |
2891 |
double *ptr_2 = &(res.getDataAtOffset(offset_2)); |
2892 |
matrix_matrix_product(SL, SM, SR, ptr_0, ptr_1, ptr_2, transpose); |
2893 |
} |
2894 |
} |
2895 |
|
2896 |
} |
2897 |
else { |
2898 |
throw DataException("Error - C_GeneralTensorProduct: unknown combination of inputs"); |
2899 |
} |
2900 |
|
2901 |
return res; |
2902 |
} |
2903 |
|
2904 |
DataAbstract* |
2905 |
Data::borrowData() const |
2906 |
{ |
2907 |
return m_data.get(); |
2908 |
} |
2909 |
|
2910 |
// Not all that happy about returning a non-const from a const |
2911 |
DataAbstract_ptr |
2912 |
Data::borrowDataPtr() const |
2913 |
{ |
2914 |
return m_data; |
2915 |
} |
2916 |
|
2917 |
// Not all that happy about returning a non-const from a const |
2918 |
DataReady_ptr |
2919 |
Data::borrowReadyPtr() const |
2920 |
{ |
2921 |
DataReady_ptr dr=dynamic_pointer_cast<DataReady>(m_data); |
2922 |
EsysAssert((dr!=0), "Error - casting to DataReady."); |
2923 |
return dr; |
2924 |
} |
2925 |
|
2926 |
std::string |
2927 |
Data::toString() const |
2928 |
{ |
2929 |
if (!m_data->isEmpty() && |
2930 |
!m_data->isLazy() && |
2931 |
getLength()>escriptParams.getInt("TOO_MANY_LINES")) |
2932 |
{ |
2933 |
stringstream temp; |
2934 |
temp << "Summary: inf="<< inf_const() << " sup=" << sup_const() << " data points=" << getNumDataPoints(); |
2935 |
return temp.str(); |
2936 |
} |
2937 |
return m_data->toString(); |
2938 |
} |
2939 |
|
2940 |
|
2941 |
|
2942 |
DataTypes::ValueType::const_reference |
2943 |
Data::getDataAtOffset(DataTypes::ValueType::size_type i) const |
2944 |
{ |
2945 |
if (isLazy()) |
2946 |
{ |
2947 |
throw DataException("Programmer error - getDataAtOffset not permitted on lazy data (object is const which prevents resolving)."); |
2948 |
} |
2949 |
return getReady()->getDataAtOffset(i); |
2950 |
} |
2951 |
|
2952 |
|
2953 |
DataTypes::ValueType::reference |
2954 |
Data::getDataAtOffset(DataTypes::ValueType::size_type i) |
2955 |
{ |
2956 |
// if (isLazy()) |
2957 |
// { |
2958 |
// throw DataException("getDataAtOffset not permitted on lazy data."); |
2959 |
// } |
2960 |
FORCERESOLVE; |
2961 |
return getReady()->getDataAtOffset(i); |
2962 |
} |
2963 |
|
2964 |
DataTypes::ValueType::const_reference |
2965 |
Data::getDataPoint(int sampleNo, int dataPointNo) const |
2966 |
{ |
2967 |
if (!isReady()) |
2968 |
{ |
2969 |
throw DataException("Programmer error - getDataPoint() not permitted on Lazy Data (object is const which prevents resolving)."); |
2970 |
} |
2971 |
else |
2972 |
{ |
2973 |
const DataReady* dr=getReady(); |
2974 |
return dr->getDataAtOffset(dr->getPointOffset(sampleNo, dataPointNo)); |
2975 |
} |
2976 |
} |
2977 |
|
2978 |
|
2979 |
DataTypes::ValueType::reference |
2980 |
Data::getDataPoint(int sampleNo, int dataPointNo) |
2981 |
{ |
2982 |
FORCERESOLVE; |
2983 |
if (!isReady()) |
2984 |
{ |
2985 |
throw DataException("Programmer error - getDataPoint() not permitted on Lazy Data."); |
2986 |
} |
2987 |
else |
2988 |
{ |
2989 |
DataReady* dr=getReady(); |
2990 |
return dr->getDataAtOffset(dr->getPointOffset(sampleNo, dataPointNo)); |
2991 |
} |
2992 |
} |
2993 |
|
2994 |
DataTypes::ValueType* |
2995 |
Data::allocSampleBuffer() const |
2996 |
{ |
2997 |
if (isLazy()) |
2998 |
{ |
2999 |
return new DataTypes::ValueType(getSampleBufferSize()); |
3000 |
} |
3001 |
else |
3002 |
{ |
3003 |
return NULL; |
3004 |
} |
3005 |
} |
3006 |
|
3007 |
void |
3008 |
Data::freeSampleBuffer(DataTypes::ValueType* buffer) |
3009 |
{ |
3010 |
if (buffer!=0) |
3011 |
{ |
3012 |
delete buffer; |
3013 |
} |
3014 |
} |
3015 |
|
3016 |
|
3017 |
/* Member functions specific to the MPI implementation */ |
3018 |
|
3019 |
void |
3020 |
Data::print() |
3021 |
{ |
3022 |
int i,j; |
3023 |
|
3024 |
printf( "Data is %dX%d\n", getNumSamples(), getNumDataPointsPerSample() ); |
3025 |
for( i=0; i<getNumSamples(); i++ ) |
3026 |
{ |
3027 |
printf( "[%6d]", i ); |
3028 |
for( j=0; j<getNumDataPointsPerSample(); j++ ) |
3029 |
printf( "\t%10.7g", (getSampleDataRW(i))[j] ); // doesn't really need RW access |
3030 |
printf( "\n" ); |
3031 |
} |
3032 |
} |
3033 |
void |
3034 |
Data::dump(const std::string fileName) const |
3035 |
{ |
3036 |
try |
3037 |
{ |
3038 |
if (isLazy()) |
3039 |
{ |
3040 |
Data temp(*this); // this is to get a non-const object which we can resolve |
3041 |
temp.resolve(); |
3042 |
temp.dump(fileName); |
3043 |
} |
3044 |
else |
3045 |
{ |
3046 |
return m_data->dump(fileName); |
3047 |
} |
3048 |
} |
3049 |
catch (exception& e) |
3050 |
{ |
3051 |
cout << e.what() << endl; |
3052 |
} |
3053 |
} |
3054 |
|
3055 |
int |
3056 |
Data::get_MPISize() const |
3057 |
{ |
3058 |
int size; |
3059 |
#ifdef PASO_MPI |
3060 |
int error; |
3061 |
error = MPI_Comm_size( get_MPIComm(), &size ); |
3062 |
#else |
3063 |
size = 1; |
3064 |
#endif |
3065 |
return size; |
3066 |
} |
3067 |
|
3068 |
int |
3069 |
Data::get_MPIRank() const |
3070 |
{ |
3071 |
int rank; |
3072 |
#ifdef PASO_MPI |
3073 |
int error; |
3074 |
error = MPI_Comm_rank( get_MPIComm(), &rank ); |
3075 |
#else |
3076 |
rank = 0; |
3077 |
#endif |
3078 |
return rank; |
3079 |
} |
3080 |
|
3081 |
MPI_Comm |
3082 |
Data::get_MPIComm() const |
3083 |
{ |
3084 |
#ifdef PASO_MPI |
3085 |
return MPI_COMM_WORLD; |
3086 |
#else |
3087 |
return -1; |
3088 |
#endif |
3089 |
} |
3090 |
|
3091 |
|