/[escript]/trunk-mpi-branch/finley/src/CPPAdapter/MeshAdapter.cpp
ViewVC logotype

Contents of /trunk-mpi-branch/finley/src/CPPAdapter/MeshAdapter.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1196 - (show annotations)
Fri Jun 15 03:45:48 2007 UTC (14 years, 5 months ago) by ksteube
File size: 56543 byte(s)
Use of PAPI on solver is now enabled with papi_instrument_solver=1 in scons/ess_options.py.
Can instrument other blocks of code with blockpapi.c.
Added interval timers to grad, integrate and Assemble_PDE.

1 // $Id$
2 /*
3 ******************************************************************************
4 * *
5 * COPYRIGHT ACcESS 2004 - All Rights Reserved *
6 * *
7 * This software is the property of ACcESS. No part of this code *
8 * may be copied in any form or by any means without the expressed written *
9 * consent of ACcESS. Copying, use or modification of this software *
10 * by any unauthorised person is illegal unless that person has a software *
11 * license agreement with ACcESS. *
12 * *
13 ******************************************************************************
14 */
15
16 #ifdef PASO_MPI
17 #include <mpi.h>
18 #endif
19 #include "MeshAdapter.h"
20
21 #include "escript/Data.h"
22 #include "escript/DataFactory.h"
23 extern "C" {
24 #include "escript/blocktimer.h"
25 }
26
27 using namespace std;
28 using namespace escript;
29
30 namespace finley {
31
32 //
33 // define the static constants
34 MeshAdapter::FunctionSpaceNamesMapType MeshAdapter::m_functionSpaceTypeNames;
35 const int MeshAdapter::DegreesOfFreedom=FINLEY_DEGREES_OF_FREEDOM;
36 const int MeshAdapter::ReducedDegreesOfFreedom=FINLEY_REDUCED_DEGREES_OF_FREEDOM;
37 const int MeshAdapter::Nodes=FINLEY_NODES;
38 const int MeshAdapter::ReducedNodes=FINLEY_REDUCED_NODES;
39 const int MeshAdapter::Elements=FINLEY_ELEMENTS;
40 const int MeshAdapter::ReducedElements=FINLEY_REDUCED_ELEMENTS;
41 const int MeshAdapter::FaceElements=FINLEY_FACE_ELEMENTS;
42 const int MeshAdapter::ReducedFaceElements=FINLEY_REDUCED_FACE_ELEMENTS;
43 const int MeshAdapter::Points=FINLEY_POINTS;
44 const int MeshAdapter::ContactElementsZero=FINLEY_CONTACT_ELEMENTS_1;
45 const int MeshAdapter::ReducedContactElementsZero=FINLEY_REDUCED_CONTACT_ELEMENTS_1;
46 const int MeshAdapter::ContactElementsOne=FINLEY_CONTACT_ELEMENTS_2;
47 const int MeshAdapter::ReducedContactElementsOne=FINLEY_REDUCED_CONTACT_ELEMENTS_2;
48
49 MeshAdapter::MeshAdapter(Finley_Mesh* finleyMesh)
50 {
51 setFunctionSpaceTypeNames();
52 //
53 // need to use a null_deleter as Finley_Mesh_dealloc deletes the pointer
54 // for us.
55 m_finleyMesh.reset(finleyMesh,null_deleter());
56 }
57
58 //
59 // The copy constructor should just increment the use count
60 MeshAdapter::MeshAdapter(const MeshAdapter& in):
61 m_finleyMesh(in.m_finleyMesh)
62 {
63 setFunctionSpaceTypeNames();
64 }
65
66 MeshAdapter::~MeshAdapter()
67 {
68 //
69 // I hope the case for the pointer being zero has been taken care of.
70 // cout << "In MeshAdapter destructor." << endl;
71 if (m_finleyMesh.unique()) {
72 // cout << "Calling dealloc." << endl;
73 Finley_Mesh_dealloc(m_finleyMesh.get());
74 // cout << "Finished dealloc." << endl;
75 }
76 }
77
78 Finley_Mesh* MeshAdapter::getFinley_Mesh() const {
79 return m_finleyMesh.get();
80 }
81
82 void MeshAdapter::write(const std::string& fileName) const
83 {
84 char *fName = (fileName.size()+1>0) ? TMPMEMALLOC(fileName.size()+1,char) : (char*)NULL;
85 strcpy(fName,fileName.c_str());
86 Finley_Mesh_write(m_finleyMesh.get(),fName);
87 checkFinleyError();
88 TMPMEMFREE(fName);
89 }
90
91 string MeshAdapter::getDescription() const
92 {
93 return "FinleyMesh";
94 }
95
96 string MeshAdapter::functionSpaceTypeAsString(int functionSpaceType) const
97 {
98 FunctionSpaceNamesMapType::iterator loc;
99 loc=m_functionSpaceTypeNames.find(functionSpaceType);
100 if (loc==m_functionSpaceTypeNames.end()) {
101 return "Invalid function space type code.";
102 } else {
103 return loc->second;
104 }
105 }
106
107 bool MeshAdapter::isValidFunctionSpaceType(int functionSpaceType) const
108 {
109 FunctionSpaceNamesMapType::iterator loc;
110 loc=m_functionSpaceTypeNames.find(functionSpaceType);
111 return (loc!=m_functionSpaceTypeNames.end());
112 }
113
114 void MeshAdapter::setFunctionSpaceTypeNames()
115 {
116 m_functionSpaceTypeNames.insert
117 (FunctionSpaceNamesMapType::value_type(DegreesOfFreedom,"Finley_DegreesOfFreedom"));
118 m_functionSpaceTypeNames.insert
119 (FunctionSpaceNamesMapType::value_type(ReducedDegreesOfFreedom,"Finley_ReducedDegreesOfFreedom"));
120 m_functionSpaceTypeNames.insert
121 (FunctionSpaceNamesMapType::value_type(Nodes,"Finley_Nodes"));
122 m_functionSpaceTypeNames.insert
123 (FunctionSpaceNamesMapType::value_type(ReducedNodes,"Finley_Reduced_Nodes"));
124 m_functionSpaceTypeNames.insert
125 (FunctionSpaceNamesMapType::value_type(Elements,"Finley_Elements"));
126 m_functionSpaceTypeNames.insert
127 (FunctionSpaceNamesMapType::value_type(ReducedElements,"Finley_Reduced_Elements"));
128 m_functionSpaceTypeNames.insert
129 (FunctionSpaceNamesMapType::value_type(FaceElements,"Finley_Face_Elements"));
130 m_functionSpaceTypeNames.insert
131 (FunctionSpaceNamesMapType::value_type(ReducedFaceElements,"Finley_Reduced_Face_Elements"));
132 m_functionSpaceTypeNames.insert
133 (FunctionSpaceNamesMapType::value_type(Points,"Finley_Points"));
134 m_functionSpaceTypeNames.insert
135 (FunctionSpaceNamesMapType::value_type(ContactElementsZero,"Finley_Contact_Elements_0"));
136 m_functionSpaceTypeNames.insert
137 (FunctionSpaceNamesMapType::value_type(ReducedContactElementsZero,"Finley_Reduced_Contact_Elements_0"));
138 m_functionSpaceTypeNames.insert
139 (FunctionSpaceNamesMapType::value_type(ContactElementsOne,"Finley_Contact_Elements_1"));
140 m_functionSpaceTypeNames.insert
141 (FunctionSpaceNamesMapType::value_type(ReducedContactElementsOne,"Finley_Reduced_Contact_Elements_1"));
142 }
143
144 int MeshAdapter::getContinuousFunctionCode() const
145 {
146 return Nodes;
147 }
148 int MeshAdapter::getReducedContinuousFunctionCode() const
149 {
150 return ReducedNodes;
151 }
152
153 int MeshAdapter::getFunctionCode() const
154 {
155 return Elements;
156 }
157 int MeshAdapter::getReducedFunctionCode() const
158 {
159 return ReducedElements;
160 }
161
162 int MeshAdapter::getFunctionOnBoundaryCode() const
163 {
164 return FaceElements;
165 }
166 int MeshAdapter::getReducedFunctionOnBoundaryCode() const
167 {
168 return ReducedFaceElements;
169 }
170
171 int MeshAdapter::getFunctionOnContactZeroCode() const
172 {
173 return ContactElementsZero;
174 }
175 int MeshAdapter::getReducedFunctionOnContactZeroCode() const
176 {
177 return ReducedContactElementsZero;
178 }
179
180 int MeshAdapter::getFunctionOnContactOneCode() const
181 {
182 return ContactElementsOne;
183 }
184 int MeshAdapter::getReducedFunctionOnContactOneCode() const
185 {
186 return ReducedContactElementsOne;
187 }
188
189 int MeshAdapter::getSolutionCode() const
190 {
191 return DegreesOfFreedom;
192 }
193
194 int MeshAdapter::getReducedSolutionCode() const
195 {
196 return ReducedDegreesOfFreedom;
197 }
198
199 int MeshAdapter::getDiracDeltaFunctionCode() const
200 {
201 return Points;
202 }
203
204 //
205 // return the spatial dimension of the Mesh:
206 //
207 int MeshAdapter::getDim() const
208 {
209 int numDim=Finley_Mesh_getDim(m_finleyMesh.get());
210 checkFinleyError();
211 return numDim;
212 }
213
214 //
215 // return the number of data points per sample and the number of samples
216 // needed to represent data on a parts of the mesh.
217 //
218 pair<int,int> MeshAdapter::getDataShape(int functionSpaceCode) const
219 {
220 int numDataPointsPerSample=0;
221 int numSamples=0;
222 Finley_Mesh* mesh=m_finleyMesh.get();
223 switch (functionSpaceCode) {
224 case(Nodes):
225 numDataPointsPerSample=1;
226 if (mesh->Nodes!=NULL) numSamples=mesh->Nodes->numNodes;
227 break;
228 case(ReducedNodes):
229 /* TODO: add ReducedNodes */
230 numDataPointsPerSample=1;
231 if (mesh->Nodes!=NULL) numSamples=mesh->Nodes->numNodes;
232 throw FinleyAdapterException("Error - ReducedNodes is not supported yet.");
233 break;
234 case(Elements):
235 if (mesh->Elements!=NULL) {
236 numSamples=mesh->Elements->numElements;
237 numDataPointsPerSample=mesh->Elements->ReferenceElement->numQuadNodes;
238 }
239 break;
240 case(ReducedElements):
241 if (mesh->Elements!=NULL) {
242 numSamples=mesh->Elements->numElements;
243 numDataPointsPerSample=mesh->Elements->ReferenceElementReducedOrder->numQuadNodes;
244 }
245 break;
246 case(FaceElements):
247 if (mesh->FaceElements!=NULL) {
248 numDataPointsPerSample=mesh->FaceElements->ReferenceElement->numQuadNodes;
249 numSamples=mesh->FaceElements->numElements;
250 }
251 break;
252 case(ReducedFaceElements):
253 if (mesh->FaceElements!=NULL) {
254 numDataPointsPerSample=mesh->FaceElements->ReferenceElementReducedOrder->numQuadNodes;
255 numSamples=mesh->FaceElements->numElements;
256 }
257 break;
258 case(Points):
259 if (mesh->Points!=NULL) {
260 numDataPointsPerSample=1;
261 numSamples=mesh->Points->numElements;
262 }
263 break;
264 case(ContactElementsZero):
265 if (mesh->ContactElements!=NULL) {
266 numDataPointsPerSample=mesh->ContactElements->ReferenceElement->numQuadNodes;
267 numSamples=mesh->ContactElements->numElements;
268 }
269 break;
270 case(ReducedContactElementsZero):
271 if (mesh->ContactElements!=NULL) {
272 numDataPointsPerSample=mesh->ContactElements->ReferenceElementReducedOrder->numQuadNodes;
273 numSamples=mesh->ContactElements->numElements;
274 }
275 break;
276 case(ContactElementsOne):
277 if (mesh->ContactElements!=NULL) {
278 numDataPointsPerSample=mesh->ContactElements->ReferenceElement->numQuadNodes;
279 numSamples=mesh->ContactElements->numElements;
280 }
281 break;
282 case(ReducedContactElementsOne):
283 if (mesh->ContactElements!=NULL) {
284 numDataPointsPerSample=mesh->ContactElements->ReferenceElementReducedOrder->numQuadNodes;
285 numSamples=mesh->ContactElements->numElements;
286 }
287 break;
288 case(DegreesOfFreedom):
289 if (mesh->Nodes!=NULL) {
290 numDataPointsPerSample=1;
291 #ifndef PASO_MPI
292 numSamples=mesh->Nodes->numDegreesOfFreedom;
293 #else
294 numSamples=mesh->Nodes->degreeOfFreedomDistribution->numLocal;
295 #endif
296 }
297 break;
298 case(ReducedDegreesOfFreedom):
299 if (mesh->Nodes!=NULL) {
300 numDataPointsPerSample=1;
301 #ifndef PASO_MPI
302 numSamples=mesh->Nodes->reducedNumDegreesOfFreedom;
303 #else
304 numSamples=mesh->Nodes->reducedDegreeOfFreedomDistribution->numLocal;
305 #endif
306 }
307 break;
308 default:
309 stringstream temp;
310 temp << "Error - Invalid function space type: " << functionSpaceCode << " for domain: " << getDescription();
311 throw FinleyAdapterException(temp.str());
312 break;
313 }
314 return pair<int,int>(numDataPointsPerSample,numSamples);
315 }
316
317 //
318 // adds linear PDE of second order into a given stiffness matrix and right hand side:
319 //
320 void MeshAdapter::addPDEToSystem(
321 SystemMatrixAdapter& mat, escript::Data& rhs,
322 const escript::Data& A, const escript::Data& B, const escript::Data& C,const escript::Data& D,const escript::Data& X,const escript::Data& Y,
323 const escript::Data& d, const escript::Data& y,
324 const escript::Data& d_contact,const escript::Data& y_contact) const
325 {
326 escriptDataC _rhs=rhs.getDataC();
327 escriptDataC _A =A.getDataC();
328 escriptDataC _B=B.getDataC();
329 escriptDataC _C=C.getDataC();
330 escriptDataC _D=D.getDataC();
331 escriptDataC _X=X.getDataC();
332 escriptDataC _Y=Y.getDataC();
333 escriptDataC _d=d.getDataC();
334 escriptDataC _y=y.getDataC();
335 escriptDataC _d_contact=d_contact.getDataC();
336 escriptDataC _y_contact=y_contact.getDataC();
337
338 Finley_Mesh* mesh=m_finleyMesh.get();
339
340 Finley_Assemble_PDE(mesh->Nodes,mesh->Elements,mat.getPaso_SystemMatrix(), &_rhs, &_A, &_B, &_C, &_D, &_X, &_Y );
341 checkFinleyError();
342
343 Finley_Assemble_PDE(mesh->Nodes,mesh->FaceElements, mat.getPaso_SystemMatrix(), &_rhs, 0, 0, 0, &_d, 0, &_y );
344 checkFinleyError();
345
346 Finley_Assemble_PDE(mesh->Nodes,mesh->ContactElements, mat.getPaso_SystemMatrix(), &_rhs , 0, 0, 0, &_d_contact, 0, &_y_contact );
347 checkFinleyError();
348 }
349
350 //
351 // adds linear PDE of second order into the right hand side only
352 //
353 void MeshAdapter::addPDEToRHS( escript::Data& rhs, const escript::Data& X,const escript::Data& Y, const escript::Data& y, const escript::Data& y_contact) const
354 {
355 Finley_Mesh* mesh=m_finleyMesh.get();
356
357 escriptDataC _rhs=rhs.getDataC();
358 escriptDataC _X=X.getDataC();
359 escriptDataC _Y=Y.getDataC();
360 escriptDataC _y=y.getDataC();
361 escriptDataC _y_contact=y_contact.getDataC();
362
363 Finley_Assemble_PDE(mesh->Nodes,mesh->Elements, 0, &_rhs, 0, 0, 0, 0, &_X, &_Y );
364 checkFinleyError();
365
366 Finley_Assemble_PDE(mesh->Nodes,mesh->FaceElements, 0, &_rhs, 0, 0, 0, 0, 0, &_y );
367 checkFinleyError();
368
369 Finley_Assemble_PDE(mesh->Nodes,mesh->ContactElements, 0, &_rhs , 0, 0, 0, 0, 0, &_y_contact );
370 checkFinleyError();
371 }
372
373 //
374 // interpolates data between different function spaces:
375 //
376 void MeshAdapter::interpolateOnDomain(escript::Data& target,const escript::Data& in) const
377 {
378 const MeshAdapter& inDomain=dynamic_cast<const MeshAdapter&>(in.getFunctionSpace().getDomain());
379 const MeshAdapter& targetDomain=dynamic_cast<const MeshAdapter&>(target.getFunctionSpace().getDomain());
380 if (inDomain!=*this)
381 throw FinleyAdapterException("Error - Illegal domain of interpolant.");
382 if (targetDomain!=*this)
383 throw FinleyAdapterException("Error - Illegal domain of interpolation target.");
384
385 Finley_Mesh* mesh=m_finleyMesh.get();
386 escriptDataC _target=target.getDataC();
387 escriptDataC _in=in.getDataC();
388 switch(in.getFunctionSpace().getTypeCode()) {
389 case(Nodes):
390 switch(target.getFunctionSpace().getTypeCode()) {
391 case(Nodes):
392 case(ReducedNodes):
393 case(DegreesOfFreedom):
394 case(ReducedDegreesOfFreedom):
395 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
396 break;
397 case(Elements):
398 case(ReducedElements):
399 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
400 break;
401 case(FaceElements):
402 case(ReducedFaceElements):
403 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
404 break;
405 case(Points):
406 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
407 break;
408 case(ContactElementsZero):
409 case(ReducedContactElementsZero):
410 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
411 break;
412 case(ContactElementsOne):
413 case(ReducedContactElementsOne):
414 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
415 break;
416 default:
417 stringstream temp;
418 temp << "Error - Interpolation on Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
419 throw FinleyAdapterException(temp.str());
420 break;
421 }
422 break;
423 case(ReducedNodes):
424 switch(target.getFunctionSpace().getTypeCode()) {
425 case(Nodes):
426 case(ReducedNodes):
427 case(DegreesOfFreedom):
428 case(ReducedDegreesOfFreedom):
429 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
430 break;
431 case(Elements):
432 case(ReducedElements):
433 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
434 break;
435 case(FaceElements):
436 case(ReducedFaceElements):
437 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
438 break;
439 case(Points):
440 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
441 break;
442 case(ContactElementsZero):
443 case(ReducedContactElementsZero):
444 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
445 break;
446 case(ContactElementsOne):
447 case(ReducedContactElementsOne):
448 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
449 break;
450 default:
451 stringstream temp;
452 temp << "Error - Interpolation on Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
453 throw FinleyAdapterException(temp.str());
454 break;
455 }
456 break;
457 case(Elements):
458 if (target.getFunctionSpace().getTypeCode()==Elements) {
459 Finley_Assemble_CopyElementData(mesh->Elements,&_target,&_in);
460 } else if (target.getFunctionSpace().getTypeCode()==ReducedElements) {
461 Finley_Assemble_AverageElementData(mesh->Elements,&_target,&_in);
462 } else {
463 throw FinleyAdapterException("Error - No interpolation with data on elements possible.");
464 }
465 break;
466 case(ReducedElements):
467 if (target.getFunctionSpace().getTypeCode()==ReducedElements) {
468 Finley_Assemble_CopyElementData(mesh->Elements,&_target,&_in);
469 } else {
470 throw FinleyAdapterException("Error - No interpolation with data on elements with reduced integration order possible.");
471 }
472 break;
473 case(FaceElements):
474 if (target.getFunctionSpace().getTypeCode()==FaceElements) {
475 Finley_Assemble_CopyElementData(mesh->FaceElements,&_target,&_in);
476 } else if (target.getFunctionSpace().getTypeCode()==ReducedFaceElements) {
477 Finley_Assemble_AverageElementData(mesh->FaceElements,&_target,&_in);
478 } else {
479 throw FinleyAdapterException("Error - No interpolation with data on face elements possible.");
480 }
481 break;
482 case(ReducedFaceElements):
483 if (target.getFunctionSpace().getTypeCode()==FaceElements) {
484 Finley_Assemble_CopyElementData(mesh->FaceElements,&_target,&_in);
485 } else {
486 throw FinleyAdapterException("Error - No interpolation with data on face elements with reduced integration order possible.");
487 }
488 break;
489 case(Points):
490 if (target.getFunctionSpace().getTypeCode()==Points) {
491 Finley_Assemble_CopyElementData(mesh->Points,&_target,&_in);
492 } else {
493 throw FinleyAdapterException("Error - No interpolation with data on points possible.");
494 }
495 break;
496 case(ContactElementsZero):
497 case(ContactElementsOne):
498 if (target.getFunctionSpace().getTypeCode()==ContactElementsZero || target.getFunctionSpace().getTypeCode()==ContactElementsOne) {
499 Finley_Assemble_CopyElementData(mesh->ContactElements,&_target,&_in);
500 } else if (target.getFunctionSpace().getTypeCode()==ReducedContactElementsZero || target.getFunctionSpace().getTypeCode()==ReducedContactElementsOne) {
501 Finley_Assemble_AverageElementData(mesh->ContactElements,&_target,&_in);
502 } else {
503 throw FinleyAdapterException("Error - No interpolation with data on contact elements possible.");
504 break;
505 }
506 break;
507 case(ReducedContactElementsZero):
508 case(ReducedContactElementsOne):
509 if (target.getFunctionSpace().getTypeCode()==ReducedContactElementsZero || target.getFunctionSpace().getTypeCode()==ReducedContactElementsOne) {
510 Finley_Assemble_CopyElementData(mesh->ContactElements,&_target,&_in);
511 } else {
512 throw FinleyAdapterException("Error - No interpolation with data on contact elements with reduced integration order possible.");
513 break;
514 }
515 break;
516 case(DegreesOfFreedom):
517 switch(target.getFunctionSpace().getTypeCode()) {
518 case(ReducedDegreesOfFreedom):
519 case(DegreesOfFreedom):
520 case(Nodes):
521 case(ReducedNodes):
522 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
523 break;
524 #ifndef PASO_MPI
525 case(Elements):
526 case(ReducedElements):
527 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
528 break;
529 case(FaceElements):
530 case(ReducedFaceElements):
531 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
532 break;
533 case(Points):
534 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
535 break;
536 case(ContactElementsZero):
537 case(ContactElementsOne):
538 case(ReducedContactElementsZero):
539 case(ReducedContactElementsOne):
540 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
541 break;
542 #else
543 /* need to copy Degrees of freedom data to nodal data so that the external values are available */
544 case(Elements):
545 case(ReducedElements):
546 { /* limit scope for variable nodeTemp */
547 escript::Data nodeTemp( in, continuousFunction(asAbstractContinuousDomain()) );
548 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&(nodeTemp.getDataC()),&_target);
549 }
550 break;
551 case(FaceElements):
552 case(ReducedFaceElements):
553 {
554 escript::Data nodeTemp( in, continuousFunction(asAbstractContinuousDomain()) );
555 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&(nodeTemp.getDataC()),&_target);
556 }
557 break;
558 case(Points):
559 {
560 escript::Data nodeTemp( in, continuousFunction(asAbstractContinuousDomain()) );
561 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&(nodeTemp.getDataC()),&_target);
562 }
563 break;
564 case(ContactElementsZero):
565 case(ContactElementsOne):
566 case(ReducedContactElementsZero):
567 case(ReducedContactElementsOne):
568 {
569 escript::Data nodeTemp( in, continuousFunction(asAbstractContinuousDomain()) );
570 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&(nodeTemp.getDataC()),&_target);
571 }
572 break;
573 #endif
574 default:
575 stringstream temp;
576 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
577 throw FinleyAdapterException(temp.str());
578 break;
579 }
580 break;
581 case(ReducedDegreesOfFreedom):
582 switch(target.getFunctionSpace().getTypeCode()) {
583 case(Nodes):
584 throw FinleyAdapterException("Error - Finley does not support interpolation from reduced degrees of freedom to mesh nodes.");
585 break;
586 case(ReducedNodes):
587 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
588 break;
589 case(DegreesOfFreedom):
590 throw FinleyAdapterException("Error - Finley does not support interpolation from reduced degrees of freedom to degrees of freedom");
591 break;
592 case(ReducedDegreesOfFreedom):
593 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
594 break;
595 case(Elements):
596 case(ReducedElements):
597 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
598 break;
599 case(FaceElements):
600 case(ReducedFaceElements):
601 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
602 break;
603 case(Points):
604 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
605 break;
606 case(ContactElementsZero):
607 case(ContactElementsOne):
608 case(ReducedContactElementsZero):
609 case(ReducedContactElementsOne):
610 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
611 break;
612 default:
613 stringstream temp;
614 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
615 throw FinleyAdapterException(temp.str());
616 break;
617 }
618 break;
619 default:
620 stringstream temp;
621 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type %d" << in.getFunctionSpace().getTypeCode();
622 throw FinleyAdapterException(temp.str());
623 break;
624 }
625 checkFinleyError();
626 }
627
628 //
629 // copies the locations of sample points into x:
630 //
631 void MeshAdapter::setToX(escript::Data& arg) const
632 {
633 const MeshAdapter& argDomain=dynamic_cast<const MeshAdapter&>(arg.getFunctionSpace().getDomain());
634 if (argDomain!=*this)
635 throw FinleyAdapterException("Error - Illegal domain of data point locations");
636 Finley_Mesh* mesh=m_finleyMesh.get();
637 // in case of values node coordinates we can do the job directly:
638 if (arg.getFunctionSpace().getTypeCode()==Nodes) {
639 escriptDataC _arg=arg.getDataC();
640 Finley_Assemble_NodeCoordinates(mesh->Nodes,&_arg);
641 } else {
642 escript::Data tmp_data=Vector(0.0,continuousFunction(asAbstractContinuousDomain()),true);
643 escriptDataC _tmp_data=tmp_data.getDataC();
644 Finley_Assemble_NodeCoordinates(mesh->Nodes,&_tmp_data);
645 // this is then interpolated onto arg:
646 interpolateOnDomain(arg,tmp_data);
647 }
648 checkFinleyError();
649 }
650
651 //
652 // return the normal vectors at the location of data points as a Data object:
653 //
654 void MeshAdapter::setToNormal(escript::Data& normal) const
655 {
656 const MeshAdapter& normalDomain=dynamic_cast<const MeshAdapter&>(normal.getFunctionSpace().getDomain());
657 if (normalDomain!=*this)
658 throw FinleyAdapterException("Error - Illegal domain of normal locations");
659 Finley_Mesh* mesh=m_finleyMesh.get();
660 escriptDataC _normal=normal.getDataC();
661 switch(normal.getFunctionSpace().getTypeCode()) {
662 case(Nodes):
663 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for nodes");
664 break;
665 case(ReducedNodes):
666 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for reduced nodes");
667 break;
668 case(Elements):
669 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for elements");
670 break;
671 case(ReducedElements):
672 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for elements with reduced integration order");
673 break;
674 case (FaceElements):
675 Finley_Assemble_setNormal(mesh->Nodes,mesh->FaceElements,&_normal);
676 break;
677 case (ReducedFaceElements):
678 Finley_Assemble_setNormal(mesh->Nodes,mesh->FaceElements,&_normal);
679 break;
680 case(Points):
681 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for point elements");
682 break;
683 case (ContactElementsOne):
684 case (ContactElementsZero):
685 Finley_Assemble_setNormal(mesh->Nodes,mesh->ContactElements,&_normal);
686 break;
687 case (ReducedContactElementsOne):
688 case (ReducedContactElementsZero):
689 Finley_Assemble_setNormal(mesh->Nodes,mesh->ContactElements,&_normal);
690 break;
691 case(DegreesOfFreedom):
692 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for degrees of freedom.");
693 break;
694 case(ReducedDegreesOfFreedom):
695 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for reduced degrees of freedom.");
696 break;
697 default:
698 stringstream temp;
699 temp << "Error - Normal Vectors: Finley does not know anything about function space type " << normal.getFunctionSpace().getTypeCode();
700 throw FinleyAdapterException(temp.str());
701 break;
702 }
703 checkFinleyError();
704 }
705
706 //
707 // interpolates data to other domain:
708 //
709 void MeshAdapter::interpolateACross(escript::Data& target,const escript::Data& source) const
710 {
711 const MeshAdapter& targetDomain=dynamic_cast<const MeshAdapter&>(target.getFunctionSpace().getDomain());
712 if (targetDomain!=*this)
713 throw FinleyAdapterException("Error - Illegal domain of interpolation target");
714
715 throw FinleyAdapterException("Error - Finley does not allow interpolation across domains yet.");
716 }
717
718 //
719 // calculates the integral of a function defined of arg:
720 //
721 void MeshAdapter::setToIntegrals(std::vector<double>& integrals,const escript::Data& arg) const
722 {
723 const MeshAdapter& argDomain=dynamic_cast<const MeshAdapter&>(arg.getFunctionSpace().getDomain());
724 if (argDomain!=*this)
725 throw FinleyAdapterException("Error - Illegal domain of integration kernel");
726
727 double blocktimer_start = blocktimer_time();
728 Finley_Mesh* mesh=m_finleyMesh.get();
729 escriptDataC _arg=arg.getDataC();
730 switch(arg.getFunctionSpace().getTypeCode()) {
731 case(Nodes):
732 /* TODO */
733 throw FinleyAdapterException("Error - Integral of data on nodes is not supported.");
734 break;
735 case(ReducedNodes):
736 /* TODO */
737 throw FinleyAdapterException("Error - Integral of data on reduced nodes is not supported.");
738 break;
739 case(Elements):
740 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_arg,&integrals[0]);
741 break;
742 case(ReducedElements):
743 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_arg,&integrals[0]);
744 break;
745 case(FaceElements):
746 Finley_Assemble_integrate(mesh->Nodes,mesh->FaceElements,&_arg,&integrals[0]);
747 break;
748 case(ReducedFaceElements):
749 Finley_Assemble_integrate(mesh->Nodes,mesh->FaceElements,&_arg,&integrals[0]);
750 break;
751 case(Points):
752 throw FinleyAdapterException("Error - Integral of data on points is not supported.");
753 break;
754 case(ContactElementsZero):
755 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
756 break;
757 case(ReducedContactElementsZero):
758 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
759 break;
760 case(ContactElementsOne):
761 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
762 break;
763 case(ReducedContactElementsOne):
764 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
765 break;
766 case(DegreesOfFreedom):
767 throw FinleyAdapterException("Error - Integral of data on degrees of freedom is not supported.");
768 break;
769 case(ReducedDegreesOfFreedom):
770 throw FinleyAdapterException("Error - Integral of data on reduced degrees of freedom is not supported.");
771 break;
772 default:
773 stringstream temp;
774 temp << "Error - Integrals: Finley does not know anything about function space type " << arg.getFunctionSpace().getTypeCode();
775 throw FinleyAdapterException(temp.str());
776 break;
777 }
778 checkFinleyError();
779 blocktimer_increment("integrate()", blocktimer_start);
780 }
781
782 //
783 // calculates the gradient of arg:
784 //
785 void MeshAdapter::setToGradient(escript::Data& grad,const escript::Data& arg) const
786 {
787 const MeshAdapter& argDomain=dynamic_cast<const MeshAdapter&>(arg.getFunctionSpace().getDomain());
788 if (argDomain!=*this)
789 throw FinleyAdapterException("Error - Illegal domain of gradient argument");
790 const MeshAdapter& gradDomain=dynamic_cast<const MeshAdapter&>(grad.getFunctionSpace().getDomain());
791 if (gradDomain!=*this)
792 throw FinleyAdapterException("Error - Illegal domain of gradient");
793
794 Finley_Mesh* mesh=m_finleyMesh.get();
795 escriptDataC _grad=grad.getDataC();
796 escriptDataC nodeDataC;
797 #ifdef PASO_MPI
798 escript::Data nodeTemp( arg, continuousFunction(asAbstractContinuousDomain()) );
799 if( arg.getFunctionSpace().getTypeCode() != Nodes )
800 {
801 Finley_Assemble_CopyNodalData(mesh->Nodes,&(nodeTemp.getDataC()),&(arg.getDataC()));
802 nodeDataC = nodeTemp.getDataC();
803 }
804 else
805 nodeDataC = arg.getDataC();
806 #else
807 nodeDataC = arg.getDataC();
808 #endif
809 switch(grad.getFunctionSpace().getTypeCode()) {
810 case(Nodes):
811 throw FinleyAdapterException("Error - Gradient at nodes is not supported.");
812 break;
813 case(ReducedNodes):
814 throw FinleyAdapterException("Error - Gradient at reduced nodes is not supported.");
815 break;
816 case(Elements):
817 Finley_Assemble_gradient(mesh->Nodes,mesh->Elements,&_grad,&nodeDataC);
818 break;
819 case(ReducedElements):
820 Finley_Assemble_gradient(mesh->Nodes,mesh->Elements,&_grad,&nodeDataC);
821 break;
822 case(FaceElements):
823 Finley_Assemble_gradient(mesh->Nodes,mesh->FaceElements,&_grad,&nodeDataC);
824 break;
825 case(ReducedFaceElements):
826 Finley_Assemble_gradient(mesh->Nodes,mesh->FaceElements,&_grad,&nodeDataC);
827 break;
828 case(Points):
829 throw FinleyAdapterException("Error - Gradient at points is not supported.");
830 break;
831 case(ContactElementsZero):
832 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
833 break;
834 case(ReducedContactElementsZero):
835 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
836 break;
837 case(ContactElementsOne):
838 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
839 break;
840 case(ReducedContactElementsOne):
841 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
842 break;
843 case(DegreesOfFreedom):
844 throw FinleyAdapterException("Error - Gradient at degrees of freedom is not supported.");
845 break;
846 case(ReducedDegreesOfFreedom):
847 throw FinleyAdapterException("Error - Gradient at reduced degrees of freedom is not supported.");
848 break;
849 default:
850 stringstream temp;
851 temp << "Error - Gradient: Finley does not know anything about function space type " << arg.getFunctionSpace().getTypeCode();
852 throw FinleyAdapterException(temp.str());
853 break;
854 }
855 checkFinleyError();
856 }
857
858 //
859 // returns the size of elements:
860 //
861 void MeshAdapter::setToSize(escript::Data& size) const
862 {
863 Finley_Mesh* mesh=m_finleyMesh.get();
864 escriptDataC tmp=size.getDataC();
865 switch(size.getFunctionSpace().getTypeCode()) {
866 case(Nodes):
867 throw FinleyAdapterException("Error - Size of nodes is not supported.");
868 break;
869 case(ReducedNodes):
870 throw FinleyAdapterException("Error - Size of reduced nodes is not supported.");
871 break;
872 case(Elements):
873 Finley_Assemble_getSize(mesh->Nodes,mesh->Elements,&tmp);
874 break;
875 case(ReducedElements):
876 Finley_Assemble_getSize(mesh->Nodes,mesh->Elements,&tmp);
877 break;
878 case(FaceElements):
879 Finley_Assemble_getSize(mesh->Nodes,mesh->FaceElements,&tmp);
880 break;
881 case(ReducedFaceElements):
882 Finley_Assemble_getSize(mesh->Nodes,mesh->FaceElements,&tmp);
883 break;
884 case(Points):
885 throw FinleyAdapterException("Error - Size of point elements is not supported.");
886 break;
887 case(ContactElementsZero):
888 case(ContactElementsOne):
889 Finley_Assemble_getSize(mesh->Nodes,mesh->ContactElements,&tmp);
890 break;
891 case(ReducedContactElementsZero):
892 case(ReducedContactElementsOne):
893 Finley_Assemble_getSize(mesh->Nodes,mesh->ContactElements,&tmp);
894 break;
895 case(DegreesOfFreedom):
896 throw FinleyAdapterException("Error - Size of degrees of freedom is not supported.");
897 break;
898 case(ReducedDegreesOfFreedom):
899 throw FinleyAdapterException("Error - Size of reduced degrees of freedom is not supported.");
900 break;
901 default:
902 stringstream temp;
903 temp << "Error - Element size: Finley does not know anything about function space type " << size.getFunctionSpace().getTypeCode();
904 throw FinleyAdapterException(temp.str());
905 break;
906 }
907 checkFinleyError();
908 }
909
910 // sets the location of nodes:
911 void MeshAdapter::setNewX(const escript::Data& new_x)
912 {
913 Finley_Mesh* mesh=m_finleyMesh.get();
914 escriptDataC tmp;
915 const MeshAdapter& newDomain=dynamic_cast<const MeshAdapter&>(new_x.getFunctionSpace().getDomain());
916 if (newDomain!=*this)
917 throw FinleyAdapterException("Error - Illegal domain of new point locations");
918 tmp = new_x.getDataC();
919 Finley_Mesh_setCoordinates(mesh,&tmp);
920 checkFinleyError();
921 }
922
923 // saves a data array in openDX format:
924 void MeshAdapter::saveDX(const std::string& filename,const boost::python::dict& arg) const
925 {
926 int MAX_namelength=256;
927 const int num_data=boost::python::extract<int>(arg.attr("__len__")());
928 /* win32 refactor */
929 char* *names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
930 for(int i=0;i<num_data;i++)
931 {
932 names[i] = (MAX_namelength>0) ? TMPMEMALLOC(MAX_namelength,char) : (char*)NULL;
933 }
934
935 char* *c_names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
936 escriptDataC *data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC) : (escriptDataC*)NULL;
937 escriptDataC* *ptr_data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC*) : (escriptDataC**)NULL;
938
939 boost::python::list keys=arg.keys();
940 for (int i=0;i<num_data;++i) {
941 escript::Data& d=boost::python::extract<escript::Data&>(arg[keys[i]]);
942 if (dynamic_cast<const MeshAdapter&>(d.getFunctionSpace().getDomain()) !=*this)
943 throw FinleyAdapterException("Error in saveDX: Data must be defined on same Domain");
944 data[i]=d.getDataC();
945 ptr_data[i]=&(data[i]);
946 std::string n=boost::python::extract<std::string>(keys[i]);
947 c_names[i]=&(names[i][0]);
948 if (n.length()>MAX_namelength-1) {
949 strncpy(c_names[i],n.c_str(),MAX_namelength-1);
950 c_names[i][MAX_namelength-1]='\0';
951 } else {
952 strcpy(c_names[i],n.c_str());
953 }
954 }
955 Finley_Mesh_saveDX(filename.c_str(),m_finleyMesh.get(),num_data,c_names,ptr_data);
956 checkFinleyError();
957
958 /* win32 refactor */
959 TMPMEMFREE(c_names);
960 TMPMEMFREE(data);
961 TMPMEMFREE(ptr_data);
962 for(int i=0;i<num_data;i++)
963 {
964 TMPMEMFREE(names[i]);
965 }
966 TMPMEMFREE(names);
967
968 return;
969 }
970
971 // saves a data array in openVTK format:
972 void MeshAdapter::saveVTK(const std::string& filename,const boost::python::dict& arg) const
973 {
974 int MAX_namelength=256;
975 const int num_data=boost::python::extract<int>(arg.attr("__len__")());
976 /* win32 refactor */
977 char* *names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
978 for(int i=0;i<num_data;i++)
979 {
980 names[i] = (MAX_namelength>0) ? TMPMEMALLOC(MAX_namelength,char) : (char*)NULL;
981 }
982
983 char* *c_names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
984 escriptDataC *data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC) : (escriptDataC*)NULL;
985 escriptDataC* *ptr_data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC*) : (escriptDataC**)NULL;
986
987 boost::python::list keys=arg.keys();
988 for (int i=0;i<num_data;++i) {
989 escript::Data& d=boost::python::extract<escript::Data&>(arg[keys[i]]);
990 if (dynamic_cast<const MeshAdapter&>(d.getFunctionSpace().getDomain()) !=*this)
991 throw FinleyAdapterException("Error in saveVTK: Data must be defined on same Domain");
992 data[i]=d.getDataC();
993 ptr_data[i]=&(data[i]);
994 std::string n=boost::python::extract<std::string>(keys[i]);
995 c_names[i]=&(names[i][0]);
996 if (n.length()>MAX_namelength-1) {
997 strncpy(c_names[i],n.c_str(),MAX_namelength-1);
998 c_names[i][MAX_namelength-1]='\0';
999 } else {
1000 strcpy(c_names[i],n.c_str());
1001 }
1002 }
1003 #ifndef PASO_MPI
1004 Finley_Mesh_saveVTK(filename.c_str(),m_finleyMesh.get(),num_data,c_names,ptr_data);
1005 #else
1006 Finley_Mesh_saveVTK_MPIO(filename.c_str(),m_finleyMesh.get(),num_data,c_names,ptr_data);
1007 #endif
1008
1009 checkFinleyError();
1010 /* win32 refactor */
1011 TMPMEMFREE(c_names);
1012 TMPMEMFREE(data);
1013 TMPMEMFREE(ptr_data);
1014 for(int i=0;i<num_data;i++)
1015 {
1016 TMPMEMFREE(names[i]);
1017 }
1018 TMPMEMFREE(names);
1019
1020 return;
1021 }
1022
1023
1024 // creates a SystemMatrixAdapter stiffness matrix an initializes it with zeros:
1025 SystemMatrixAdapter MeshAdapter::newSystemMatrix(
1026 const int row_blocksize,
1027 const escript::FunctionSpace& row_functionspace,
1028 const int column_blocksize,
1029 const escript::FunctionSpace& column_functionspace,
1030 const int type) const
1031 {
1032 int reduceRowOrder=0;
1033 int reduceColOrder=0;
1034 // is the domain right?
1035 const MeshAdapter& row_domain=dynamic_cast<const MeshAdapter&>(row_functionspace.getDomain());
1036 if (row_domain!=*this)
1037 throw FinleyAdapterException("Error - domain of row function space does not match the domain of matrix generator.");
1038 const MeshAdapter& col_domain=dynamic_cast<const MeshAdapter&>(column_functionspace.getDomain());
1039 if (col_domain!=*this)
1040 throw FinleyAdapterException("Error - domain of columnn function space does not match the domain of matrix generator.");
1041 // is the function space type right
1042 if (row_functionspace.getTypeCode()==DegreesOfFreedom) {
1043 reduceRowOrder=0;
1044 } else if (row_functionspace.getTypeCode()==ReducedDegreesOfFreedom) {
1045 reduceRowOrder=1;
1046 } else {
1047 throw FinleyAdapterException("Error - illegal function space type for system matrix rows.");
1048 }
1049 if (column_functionspace.getTypeCode()==DegreesOfFreedom) {
1050 reduceColOrder=0;
1051 } else if (column_functionspace.getTypeCode()==ReducedDegreesOfFreedom) {
1052 reduceColOrder=1;
1053 } else {
1054 throw FinleyAdapterException("Error - illegal function space type for system matrix columns.");
1055 }
1056 // generate matrix:
1057
1058 Paso_SystemMatrixPattern* fsystemMatrixPattern=Finley_getPattern(getFinley_Mesh(),reduceRowOrder,reduceColOrder);
1059 checkFinleyError();
1060 Paso_SystemMatrix* fsystemMatrix;
1061 int trilinos = 0;
1062 if (trilinos) {
1063 #ifdef TRILINOS
1064 /* Allocation Epetra_VrbMatrix here */
1065 #endif
1066 }
1067 else {
1068 fsystemMatrix=Paso_SystemMatrix_alloc(type,fsystemMatrixPattern,row_blocksize,column_blocksize);
1069 }
1070 checkPasoError();
1071 Paso_SystemMatrixPattern_dealloc(fsystemMatrixPattern);
1072 return SystemMatrixAdapter(fsystemMatrix,row_blocksize,row_functionspace,column_blocksize,column_functionspace);
1073 }
1074
1075 //
1076 // vtkObject MeshAdapter::createVtkObject() const
1077 // TODO:
1078 //
1079
1080 //
1081 // returns true if data at the atom_type is considered as being cell centered:
1082 bool MeshAdapter::isCellOriented(int functionSpaceCode) const
1083 {
1084 switch(functionSpaceCode) {
1085 case(Nodes):
1086 case(DegreesOfFreedom):
1087 case(ReducedDegreesOfFreedom):
1088 return false;
1089 break;
1090 case(Elements):
1091 case(FaceElements):
1092 case(Points):
1093 case(ContactElementsZero):
1094 case(ContactElementsOne):
1095 case(ReducedElements):
1096 case(ReducedFaceElements):
1097 case(ReducedContactElementsZero):
1098 case(ReducedContactElementsOne):
1099 return true;
1100 break;
1101 default:
1102 stringstream temp;
1103 temp << "Error - Cell: Finley does not know anything about function space type " << functionSpaceCode;
1104 throw FinleyAdapterException(temp.str());
1105 break;
1106 }
1107 checkFinleyError();
1108 return false;
1109 }
1110
1111 bool MeshAdapter::probeInterpolationOnDomain(int functionSpaceType_source,int functionSpaceType_target) const
1112 {
1113 switch(functionSpaceType_source) {
1114 case(Nodes):
1115 switch(functionSpaceType_target) {
1116 case(Nodes):
1117 case(ReducedNodes):
1118 case(ReducedDegreesOfFreedom):
1119 case(DegreesOfFreedom):
1120 case(Elements):
1121 case(ReducedElements):
1122 case(FaceElements):
1123 case(ReducedFaceElements):
1124 case(Points):
1125 case(ContactElementsZero):
1126 case(ReducedContactElementsZero):
1127 case(ContactElementsOne):
1128 case(ReducedContactElementsOne):
1129 return true;
1130 default:
1131 stringstream temp;
1132 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1133 throw FinleyAdapterException(temp.str());
1134 }
1135 break;
1136 case(ReducedNodes):
1137 switch(functionSpaceType_target) {
1138 case(ReducedNodes):
1139 case(ReducedDegreesOfFreedom):
1140 case(Elements):
1141 case(ReducedElements):
1142 case(FaceElements):
1143 case(ReducedFaceElements):
1144 case(Points):
1145 case(ContactElementsZero):
1146 case(ReducedContactElementsZero):
1147 case(ContactElementsOne):
1148 case(ReducedContactElementsOne):
1149 return true;
1150 case(Nodes):
1151 case(DegreesOfFreedom):
1152 return false;
1153 default:
1154 stringstream temp;
1155 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1156 throw FinleyAdapterException(temp.str());
1157 }
1158 break;
1159 case(Elements):
1160 if (functionSpaceType_target==Elements) {
1161 return true;
1162 } else if (functionSpaceType_target==ReducedElements) {
1163 return true;
1164 } else {
1165 return false;
1166 }
1167 case(ReducedElements):
1168 if (functionSpaceType_target==ReducedElements) {
1169 return true;
1170 } else {
1171 return false;
1172 }
1173 case(FaceElements):
1174 if (functionSpaceType_target==FaceElements) {
1175 return true;
1176 } else if (functionSpaceType_target==ReducedFaceElements) {
1177 return true;
1178 } else {
1179 return false;
1180 }
1181 case(ReducedFaceElements):
1182 if (functionSpaceType_target==ReducedFaceElements) {
1183 return true;
1184 } else {
1185 return false;
1186 }
1187 case(Points):
1188 if (functionSpaceType_target==Points) {
1189 return true;
1190 } else {
1191 return false;
1192 }
1193 case(ContactElementsZero):
1194 case(ContactElementsOne):
1195 if (functionSpaceType_target==ContactElementsZero || functionSpaceType_target==ContactElementsOne) {
1196 return true;
1197 } else if (functionSpaceType_target==ReducedContactElementsZero || functionSpaceType_target==ReducedContactElementsOne) {
1198 return true;
1199 } else {
1200 return false;
1201 }
1202 case(ReducedContactElementsZero):
1203 case(ReducedContactElementsOne):
1204 if (functionSpaceType_target==ReducedContactElementsZero || functionSpaceType_target==ReducedContactElementsOne) {
1205 return true;
1206 } else {
1207 return false;
1208 }
1209 case(DegreesOfFreedom):
1210 switch(functionSpaceType_target) {
1211 case(ReducedDegreesOfFreedom):
1212 case(DegreesOfFreedom):
1213 case(Nodes):
1214 case(ReducedNodes):
1215 case(Elements):
1216 case(ReducedElements):
1217 case(Points):
1218 case(FaceElements):
1219 case(ReducedFaceElements):
1220 case(ContactElementsZero):
1221 case(ReducedContactElementsZero):
1222 case(ContactElementsOne):
1223 case(ReducedContactElementsOne):
1224 return true;
1225 default:
1226 stringstream temp;
1227 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1228 throw FinleyAdapterException(temp.str());
1229 }
1230 break;
1231 case(ReducedDegreesOfFreedom):
1232 switch(functionSpaceType_target) {
1233 case(ReducedDegreesOfFreedom):
1234 case(ReducedNodes):
1235 case(Elements):
1236 case(ReducedElements):
1237 case(FaceElements):
1238 case(ReducedFaceElements):
1239 case(Points):
1240 case(ContactElementsZero):
1241 case(ReducedContactElementsZero):
1242 case(ContactElementsOne):
1243 case(ReducedContactElementsOne):
1244 return true;
1245 case(Nodes):
1246 case(DegreesOfFreedom):
1247 return false;
1248 default:
1249 stringstream temp;
1250 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1251 throw FinleyAdapterException(temp.str());
1252 }
1253 break;
1254 default:
1255 stringstream temp;
1256 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_source;
1257 throw FinleyAdapterException(temp.str());
1258 break;
1259 }
1260 checkFinleyError();
1261 return false;
1262 }
1263
1264 bool MeshAdapter::probeInterpolationACross(int functionSpaceType_source,const AbstractDomain& targetDomain, int functionSpaceType_target) const
1265 {
1266 return false;
1267 }
1268
1269 bool MeshAdapter::operator==(const AbstractDomain& other) const
1270 {
1271 const MeshAdapter* temp=dynamic_cast<const MeshAdapter*>(&other);
1272 if (temp!=0) {
1273 return (m_finleyMesh==temp->m_finleyMesh);
1274 } else {
1275 return false;
1276 }
1277 }
1278
1279 bool MeshAdapter::operator!=(const AbstractDomain& other) const
1280 {
1281 return !(operator==(other));
1282 }
1283
1284 int MeshAdapter::getSystemMatrixTypeId(const int solver, const int package, const bool symmetry) const
1285 {
1286 int out=Paso_SystemMatrix_getSystemMatrixTypeId(SystemMatrixAdapter::mapOptionToPaso(solver),SystemMatrixAdapter::mapOptionToPaso(package),symmetry?1:0);
1287 checkPasoError();
1288 return out;
1289 }
1290
1291 escript::Data MeshAdapter::getX() const
1292 {
1293 return continuousFunction(asAbstractContinuousDomain()).getX();
1294 }
1295
1296 escript::Data MeshAdapter::getNormal() const
1297 {
1298 return functionOnBoundary(asAbstractContinuousDomain()).getNormal();
1299 }
1300
1301 escript::Data MeshAdapter::getSize() const
1302 {
1303 return function(asAbstractContinuousDomain()).getSize();
1304 }
1305
1306 int* MeshAdapter::borrowSampleReferenceIDs(int functionSpaceType) const
1307 {
1308 int *out=0,i;
1309 Finley_Mesh* mesh=m_finleyMesh.get();
1310 switch (functionSpaceType) {
1311 case(Nodes):
1312 if (mesh->Nodes!=NULL) {
1313 out=mesh->Nodes->Id;
1314 break;
1315 }
1316 case(ReducedNodes):
1317 throw FinleyAdapterException("Error - ReducedNodes not supported yet.");
1318 break;
1319 case(Elements):
1320 out=mesh->Elements->Id;
1321 break;
1322 case(ReducedElements):
1323 out=mesh->Elements->Id;
1324 break;
1325 case(FaceElements):
1326 out=mesh->FaceElements->Id;
1327 break;
1328 case(ReducedFaceElements):
1329 out=mesh->FaceElements->Id;
1330 break;
1331 case(Points):
1332 out=mesh->Points->Id;
1333 break;
1334 case(ContactElementsZero):
1335 out=mesh->ContactElements->Id;
1336 break;
1337 case(ReducedContactElementsZero):
1338 out=mesh->ContactElements->Id;
1339 break;
1340 case(ContactElementsOne):
1341 out=mesh->ContactElements->Id;
1342 break;
1343 case(ReducedContactElementsOne):
1344 out=mesh->ContactElements->Id;
1345 break;
1346 case(DegreesOfFreedom):
1347 out=mesh->Nodes->degreeOfFreedomId;
1348 break;
1349 case(ReducedDegreesOfFreedom):
1350 out=mesh->Nodes->reducedDegreeOfFreedomId;
1351 break;
1352 default:
1353 stringstream temp;
1354 temp << "Error - Invalid function space type: " << functionSpaceType << " for domain: " << getDescription();
1355 throw FinleyAdapterException(temp.str());
1356 break;
1357 }
1358 return out;
1359 }
1360 int MeshAdapter::getTagFromSampleNo(int functionSpaceType, int sampleNo) const
1361 {
1362 int out=0;
1363 Finley_Mesh* mesh=m_finleyMesh.get();
1364 switch (functionSpaceType) {
1365 case(Nodes):
1366 out=mesh->Nodes->Tag[sampleNo];
1367 break;
1368 case(ReducedNodes):
1369 throw FinleyAdapterException(" Error - ReducedNodes does not support tags.");
1370 break;
1371 case(Elements):
1372 out=mesh->Elements->Tag[sampleNo];
1373 break;
1374 case(ReducedElements):
1375 out=mesh->Elements->Tag[sampleNo];
1376 break;
1377 case(FaceElements):
1378 out=mesh->FaceElements->Tag[sampleNo];
1379 break;
1380 case(ReducedFaceElements):
1381 out=mesh->FaceElements->Tag[sampleNo];
1382 break;
1383 case(Points):
1384 out=mesh->Points->Tag[sampleNo];
1385 break;
1386 case(ContactElementsZero):
1387 out=mesh->ContactElements->Tag[sampleNo];
1388 break;
1389 case(ReducedContactElementsZero):
1390 out=mesh->ContactElements->Tag[sampleNo];
1391 break;
1392 case(ContactElementsOne):
1393 out=mesh->ContactElements->Tag[sampleNo];
1394 break;
1395 case(ReducedContactElementsOne):
1396 out=mesh->ContactElements->Tag[sampleNo];
1397 break;
1398 case(DegreesOfFreedom):
1399 throw FinleyAdapterException(" Error - DegreesOfFreedom does not support tags.");
1400 break;
1401 case(ReducedDegreesOfFreedom):
1402 throw FinleyAdapterException(" Error - ReducedDegreesOfFreedom does not support tags.");
1403 break;
1404 default:
1405 stringstream temp;
1406 temp << "Error - Invalid function space type: " << functionSpaceType << " for domain: " << getDescription();
1407 throw FinleyAdapterException(temp.str());
1408 break;
1409 }
1410 return out;
1411 }
1412
1413
1414 void MeshAdapter::setTags(const int functionSpaceType, const int newTag, const escript::Data& mask) const
1415 {
1416 Finley_Mesh* mesh=m_finleyMesh.get();
1417 escriptDataC tmp=mask.getDataC();
1418 switch(functionSpaceType) {
1419 case(Nodes):
1420 Finley_NodeFile_setTags(mesh->Nodes,newTag,&tmp);
1421 break;
1422 case(ReducedNodes):
1423 throw FinleyAdapterException("Error - ReducedNodes does not support tags");
1424 break;
1425 case(DegreesOfFreedom):
1426 throw FinleyAdapterException("Error - DegreesOfFreedom does not support tags");
1427 break;
1428 case(ReducedDegreesOfFreedom):
1429 throw FinleyAdapterException("Error - ReducedDegreesOfFreedom does not support tags");
1430 break;
1431 case(Elements):
1432 Finley_ElementFile_setTags(mesh->Elements,newTag,&tmp);
1433 break;
1434 case(ReducedElements):
1435 Finley_ElementFile_setTags(mesh->Elements,newTag,&tmp);
1436 break;
1437 case(FaceElements):
1438 Finley_ElementFile_setTags(mesh->FaceElements,newTag,&tmp);
1439 break;
1440 case(ReducedFaceElements):
1441 Finley_ElementFile_setTags(mesh->FaceElements,newTag,&tmp);
1442 break;
1443 case(Points):
1444 Finley_ElementFile_setTags(mesh->Points,newTag,&tmp);
1445 break;
1446 case(ContactElementsZero):
1447 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
1448 break;
1449 case(ReducedContactElementsZero):
1450 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
1451 break;
1452 case(ContactElementsOne):
1453 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
1454 break;
1455 case(ReducedContactElementsOne):
1456 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
1457 break;
1458 default:
1459 stringstream temp;
1460 temp << "Error - Finley does not know anything about function space type " << functionSpaceType;
1461 throw FinleyAdapterException(temp.str());
1462 }
1463 checkFinleyError();
1464 return;
1465 }
1466
1467 void MeshAdapter::setTagMap(const std::string& name, int tag)
1468 {
1469 Finley_Mesh* mesh=m_finleyMesh.get();
1470 Finley_Mesh_addTagMap(mesh, name.c_str(),tag);
1471 checkPasoError();
1472 // throwStandardException("MeshAdapter::set TagMap is not implemented.");
1473 }
1474
1475 int MeshAdapter::getTag(const std::string& name) const
1476 {
1477 Finley_Mesh* mesh=m_finleyMesh.get();
1478 int tag=0;
1479 tag=Finley_Mesh_getTag(mesh, name.c_str());
1480 checkPasoError();
1481 // throwStandardException("MeshAdapter::getTag is not implemented.");
1482 return tag;
1483 }
1484
1485 bool MeshAdapter::isValidTagName(const std::string& name) const
1486 {
1487 Finley_Mesh* mesh=m_finleyMesh.get();
1488 return Finley_Mesh_isValidTagName(mesh,name.c_str());
1489 }
1490
1491 std::string MeshAdapter::showTagNames() const
1492 {
1493 stringstream temp;
1494 Finley_Mesh* mesh=m_finleyMesh.get();
1495 Finley_TagMap* tag_map=mesh->TagMap;
1496 while (tag_map) {
1497 temp << tag_map->name;
1498 tag_map=tag_map->next;
1499 if (tag_map) temp << ", ";
1500 }
1501 return temp.str();
1502 }
1503
1504 } // end of namespace

Properties

Name Value
svn:eol-style native
svn:keywords Author Date Id Revision

  ViewVC Help
Powered by ViewVC 1.1.26