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

Annotation of /trunk/finley/src/CPPAdapter/MeshAdapter.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1064 - (hide annotations)
Tue Mar 27 06:21:02 2007 UTC (16 years ago) by gross
File size: 55232 byte(s)
test for reduced integration order for grad, interpolate and integrate added.
Bug shown by the tests have been fixed.


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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26