1 |
|
2 |
/******************************************************* |
3 |
* |
4 |
* Copyright (c) 2003-2009 by University of Queensland |
5 |
* Earth Systems Science Computational Center (ESSCC) |
6 |
* http://www.uq.edu.au/esscc |
7 |
* |
8 |
* Primary Business: Queensland, Australia |
9 |
* Licensed under the Open Software License version 3.0 |
10 |
* http://www.opensource.org/licenses/osl-3.0.php |
11 |
* |
12 |
*******************************************************/ |
13 |
|
14 |
|
15 |
#if !defined finley_MeshAdapter_20040526_H |
16 |
#define finley_MeshAdapter_20040526_H |
17 |
#include "system_dep.h" |
18 |
|
19 |
extern "C" { |
20 |
#include "../Mesh.h" |
21 |
#include "../Finley.h" |
22 |
#include "../Assemble.h" |
23 |
#include "paso/SystemMatrix.h" |
24 |
#include "paso/SolverFCT.h" |
25 |
#include "paso/Paso_MPI.h" |
26 |
} |
27 |
|
28 |
#include "FinleyError.h" |
29 |
#include "FinleyAdapterException.h" |
30 |
|
31 |
#include "SystemMatrixAdapter.h" |
32 |
#include "TransportProblemAdapter.h" |
33 |
#include "escript/AbstractContinuousDomain.h" |
34 |
#include "escript/FunctionSpace.h" |
35 |
#include "escript/FunctionSpaceFactory.h" |
36 |
|
37 |
#include <boost/shared_ptr.hpp> |
38 |
#include <boost/python/dict.hpp> |
39 |
#include <boost/python/extract.hpp> |
40 |
|
41 |
#include <map> |
42 |
#include <vector> |
43 |
#include <string> |
44 |
#include <sstream> |
45 |
|
46 |
// |
47 |
// forward declarations |
48 |
class Data; |
49 |
|
50 |
//using namespace escript; |
51 |
|
52 |
namespace finley { |
53 |
|
54 |
struct null_deleter |
55 |
{ |
56 |
void operator()(void const *ptr) const |
57 |
{ |
58 |
} |
59 |
}; |
60 |
|
61 |
|
62 |
/** |
63 |
\brief |
64 |
MeshAdapter implements the AbstractContinuousDomain |
65 |
interface for the Finley library. |
66 |
|
67 |
Description: |
68 |
MeshAdapter implements the AbstractContinuousDomain |
69 |
interface for the Finley library. |
70 |
*/ |
71 |
|
72 |
class MeshAdapter : public escript::AbstractContinuousDomain { |
73 |
|
74 |
public: |
75 |
|
76 |
// |
77 |
// Codes for function space types supported |
78 |
static const int DegreesOfFreedom; |
79 |
static const int ReducedDegreesOfFreedom; |
80 |
static const int Nodes; |
81 |
static const int ReducedNodes; |
82 |
static const int Elements; |
83 |
static const int ReducedElements; |
84 |
static const int FaceElements; |
85 |
static const int ReducedFaceElements; |
86 |
static const int Points; |
87 |
static const int ContactElementsZero; |
88 |
static const int ReducedContactElementsZero; |
89 |
static const int ContactElementsOne; |
90 |
static const int ReducedContactElementsOne; |
91 |
|
92 |
/** |
93 |
\brief |
94 |
Constructor for MeshAdapter |
95 |
|
96 |
Description: |
97 |
Constructor for MeshAdapter. The pointer passed to MeshAdapter |
98 |
is deleted using a call to Finley_Mesh_free in the |
99 |
MeshAdapter destructor. |
100 |
|
101 |
Throws: |
102 |
May throw an exception derived from EsysException |
103 |
|
104 |
\param finleyMesh Input - A pointer to the externally constructed |
105 |
finley mesh.The pointer passed to MeshAdapter |
106 |
is deleted using a call to |
107 |
Finley_Mesh_free in the MeshAdapter |
108 |
destructor. |
109 |
*/ |
110 |
FINLEY_DLL_API |
111 |
MeshAdapter(Finley_Mesh* finleyMesh=0); |
112 |
|
113 |
/** |
114 |
\brief |
115 |
Copy constructor. |
116 |
*/ |
117 |
FINLEY_DLL_API |
118 |
MeshAdapter(const MeshAdapter& in); |
119 |
|
120 |
/** |
121 |
\brief |
122 |
Destructor for MeshAdapter. As specified in the constructor |
123 |
this calls Finley_Mesh_free for the pointer given to the |
124 |
constructor. |
125 |
*/ |
126 |
FINLEY_DLL_API |
127 |
~MeshAdapter(); |
128 |
|
129 |
/** |
130 |
\brief |
131 |
return the number of processors used for this domain |
132 |
*/ |
133 |
FINLEY_DLL_API |
134 |
virtual int getMPISize() const; |
135 |
/** |
136 |
\brief |
137 |
return the number MPI rank of this processor |
138 |
*/ |
139 |
|
140 |
FINLEY_DLL_API |
141 |
virtual int getMPIRank() const; |
142 |
|
143 |
/** |
144 |
\brief |
145 |
If compiled for MPI then execute an MPI_Barrier, else do nothing |
146 |
*/ |
147 |
|
148 |
FINLEY_DLL_API |
149 |
virtual void MPIBarrier() const; |
150 |
|
151 |
/** |
152 |
\brief |
153 |
Return true if on MPI processor 0, else false |
154 |
*/ |
155 |
|
156 |
FINLEY_DLL_API |
157 |
virtual bool onMasterProcessor() const; |
158 |
|
159 |
FINLEY_DLL_API |
160 |
#ifdef PASO_MPI |
161 |
MPI_Comm |
162 |
#else |
163 |
unsigned int |
164 |
#endif |
165 |
getMPIComm() const; |
166 |
|
167 |
/** |
168 |
\brief |
169 |
return this as an AbstractContinuousDomain. |
170 |
*/ |
171 |
inline const AbstractContinuousDomain& asAbstractContinuousDomain() const |
172 |
{ |
173 |
return *(static_cast<const AbstractContinuousDomain*>(this)); |
174 |
} |
175 |
|
176 |
/** |
177 |
\brief |
178 |
return this as an AbstractDomain. |
179 |
*/ |
180 |
inline const AbstractDomain& asAbstractDomain() const |
181 |
{ |
182 |
return *(static_cast<const AbstractDomain*>(this)); |
183 |
} |
184 |
|
185 |
/** |
186 |
\brief |
187 |
Write the current mesh to a file with the given name. |
188 |
\param fileName Input - The name of the file to write to. |
189 |
*/ |
190 |
FINLEY_DLL_API |
191 |
void write(const std::string& fileName) const; |
192 |
|
193 |
/** |
194 |
\brief |
195 |
\param full |
196 |
*/ |
197 |
FINLEY_DLL_API |
198 |
void Print_Mesh_Info(const bool full=false) const; |
199 |
|
200 |
/** |
201 |
\brief |
202 |
dumps the mesh to a file with the given name. |
203 |
\param fileName Input - The name of the file |
204 |
*/ |
205 |
FINLEY_DLL_API |
206 |
void dump(const std::string& fileName) const; |
207 |
|
208 |
/** |
209 |
\brief |
210 |
return the pointer to the underlying finley mesh structure |
211 |
*/ |
212 |
FINLEY_DLL_API |
213 |
Finley_Mesh* getFinley_Mesh() const; |
214 |
|
215 |
/** |
216 |
\brief |
217 |
Return the tag key for the given sample number. |
218 |
\param functionSpaceType Input - The function space type. |
219 |
\param sampleNo Input - The sample number. |
220 |
*/ |
221 |
FINLEY_DLL_API |
222 |
int getTagFromSampleNo(int functionSpaceType, int sampleNo) const; |
223 |
|
224 |
/** |
225 |
\brief |
226 |
Return the reference number of the given sample number. |
227 |
\param functionSpaceType Input - The function space type. |
228 |
*/ |
229 |
FINLEY_DLL_API |
230 |
const int* borrowSampleReferenceIDs(int functionSpaceType) const; |
231 |
|
232 |
/** |
233 |
\brief |
234 |
Returns true if the given integer is a valid function space type |
235 |
for this domain. |
236 |
*/ |
237 |
FINLEY_DLL_API |
238 |
virtual bool isValidFunctionSpaceType(int functionSpaceType) const; |
239 |
|
240 |
/** |
241 |
\brief |
242 |
Return a description for this domain |
243 |
*/ |
244 |
FINLEY_DLL_API |
245 |
virtual std::string getDescription() const; |
246 |
|
247 |
/** |
248 |
\brief |
249 |
Return a description for the given function space type code |
250 |
*/ |
251 |
FINLEY_DLL_API |
252 |
virtual std::string functionSpaceTypeAsString(int functionSpaceType) const; |
253 |
|
254 |
/** |
255 |
\brief |
256 |
Build the table of function space type names |
257 |
*/ |
258 |
FINLEY_DLL_API |
259 |
void setFunctionSpaceTypeNames(); |
260 |
|
261 |
/** |
262 |
\brief |
263 |
Return a continuous FunctionSpace code |
264 |
*/ |
265 |
FINLEY_DLL_API |
266 |
virtual int getContinuousFunctionCode() const; |
267 |
|
268 |
/** |
269 |
\brief |
270 |
Return a continuous on reduced order nodes FunctionSpace code |
271 |
*/ |
272 |
FINLEY_DLL_API |
273 |
virtual int getReducedContinuousFunctionCode() const; |
274 |
|
275 |
/** |
276 |
\brief |
277 |
Return a function FunctionSpace code |
278 |
*/ |
279 |
FINLEY_DLL_API |
280 |
virtual int getFunctionCode() const; |
281 |
|
282 |
/** |
283 |
\brief |
284 |
Return a function with reduced integration order FunctionSpace code |
285 |
*/ |
286 |
FINLEY_DLL_API |
287 |
virtual int getReducedFunctionCode() const; |
288 |
|
289 |
/** |
290 |
\brief |
291 |
Return a function on boundary FunctionSpace code |
292 |
*/ |
293 |
FINLEY_DLL_API |
294 |
virtual int getFunctionOnBoundaryCode() const; |
295 |
|
296 |
/** |
297 |
\brief |
298 |
Return a function on boundary with reduced integration order FunctionSpace code |
299 |
*/ |
300 |
FINLEY_DLL_API |
301 |
virtual int getReducedFunctionOnBoundaryCode() const; |
302 |
|
303 |
/** |
304 |
\brief |
305 |
Return a FunctionOnContactZero code |
306 |
*/ |
307 |
FINLEY_DLL_API |
308 |
virtual int getFunctionOnContactZeroCode() const; |
309 |
|
310 |
/** |
311 |
\brief |
312 |
Return a FunctionOnContactZero code with reduced integration order |
313 |
*/ |
314 |
FINLEY_DLL_API |
315 |
virtual int getReducedFunctionOnContactZeroCode() const; |
316 |
|
317 |
/** |
318 |
\brief |
319 |
Return a FunctionOnContactOne code |
320 |
*/ |
321 |
FINLEY_DLL_API |
322 |
virtual int getFunctionOnContactOneCode() const; |
323 |
|
324 |
/** |
325 |
\brief |
326 |
Return a FunctionOnContactOne code with reduced integration order |
327 |
*/ |
328 |
FINLEY_DLL_API |
329 |
virtual int getReducedFunctionOnContactOneCode() const; |
330 |
|
331 |
/** |
332 |
\brief |
333 |
Return a Solution code |
334 |
*/ |
335 |
FINLEY_DLL_API |
336 |
virtual int getSolutionCode() const; |
337 |
|
338 |
/** |
339 |
\brief |
340 |
Return a ReducedSolution code |
341 |
*/ |
342 |
FINLEY_DLL_API |
343 |
virtual int getReducedSolutionCode() const; |
344 |
|
345 |
/** |
346 |
\brief |
347 |
Return a DiracDeltaFunction code |
348 |
*/ |
349 |
FINLEY_DLL_API |
350 |
virtual int getDiracDeltaFunctionCode() const; |
351 |
|
352 |
/** |
353 |
5B |
354 |
\brief |
355 |
*/ |
356 |
typedef std::map<int, std::string> FunctionSpaceNamesMapType; |
357 |
|
358 |
/** |
359 |
\brief |
360 |
*/ |
361 |
FINLEY_DLL_API |
362 |
virtual int getDim() const; |
363 |
|
364 |
/** |
365 |
\brief |
366 |
Returns a status indicator of the domain. The status identifier should be unique over |
367 |
the live time if the object but may be updated if changes to the domain happen, e.g. |
368 |
modifications to its geometry. |
369 |
|
370 |
This has to be implemented by the actual Domain adapter. |
371 |
*/ |
372 |
FINLEY_DLL_API |
373 |
virtual StatusType getStatus() const; |
374 |
|
375 |
|
376 |
/** |
377 |
\brief |
378 |
Return the number of data points summed across all MPI processes |
379 |
*/ |
380 |
FINLEY_DLL_API |
381 |
virtual int getNumDataPointsGlobal() const; |
382 |
|
383 |
/** |
384 |
\brief |
385 |
Return the number of data points per sample, and the number of samples as a pair. |
386 |
\param functionSpaceCode Input - |
387 |
*/ |
388 |
FINLEY_DLL_API |
389 |
virtual std::pair<int,int> getDataShape(int functionSpaceCode) const; |
390 |
|
391 |
/** |
392 |
\brief |
393 |
copies the location of data points into arg. The domain of arg has to match this. |
394 |
has to be implemented by the actual Domain adapter. |
395 |
*/ |
396 |
FINLEY_DLL_API |
397 |
virtual void setToX(escript::Data& arg) const; |
398 |
|
399 |
/** |
400 |
\brief |
401 |
sets a map from a clear tag name to a tag key |
402 |
\param name Input - tag name. |
403 |
\param tag Input - tag key. |
404 |
*/ |
405 |
FINLEY_DLL_API |
406 |
virtual void setTagMap(const std::string& name, int tag); |
407 |
|
408 |
/** |
409 |
\brief |
410 |
Return the tag key for tag name. |
411 |
\param name Input - tag name |
412 |
*/ |
413 |
FINLEY_DLL_API |
414 |
virtual int getTag(const std::string& name) const; |
415 |
|
416 |
/** |
417 |
\brief |
418 |
Returns true if name is a defined tage name. |
419 |
\param name Input - tag name to be checked. |
420 |
*/ |
421 |
FINLEY_DLL_API |
422 |
virtual bool isValidTagName(const std::string& name) const; |
423 |
|
424 |
/** |
425 |
\brief |
426 |
Returns all tag names in a single string sperated by commas |
427 |
*/ |
428 |
FINLEY_DLL_API |
429 |
virtual std::string showTagNames() const; |
430 |
|
431 |
/** |
432 |
\brief |
433 |
assigns new location to the domain |
434 |
*/ |
435 |
FINLEY_DLL_API |
436 |
virtual void setNewX(const escript::Data& arg); |
437 |
|
438 |
/** |
439 |
\brief |
440 |
interpolates data given on source onto target where source and target have to be given on the same domain. |
441 |
*/ |
442 |
FINLEY_DLL_API |
443 |
virtual void interpolateOnDomain(escript::Data& target,const escript::Data& source) const; |
444 |
|
445 |
|
446 |
FINLEY_DLL_API |
447 |
virtual bool probeInterpolationOnDomain(int functionSpaceType_source,int functionSpaceType_target) const; |
448 |
|
449 |
/** |
450 |
\brief given a vector of FunctionSpace typecodes, pass back a code which then can all be interpolated to. |
451 |
\return true is result is valid, false if not |
452 |
*/ |
453 |
FINLEY_DLL_API |
454 |
bool |
455 |
commonFunctionSpace(const std::vector<int>& fs, int& resultcode) const; |
456 |
|
457 |
/** |
458 |
\brief |
459 |
interpolates data given on source onto target where source and target are given on different domains. |
460 |
has to be implemented by the actual Domain adapter. |
461 |
*/ |
462 |
FINLEY_DLL_API |
463 |
virtual void interpolateACross(escript::Data& target, const escript::Data& source) const; |
464 |
|
465 |
/** |
466 |
\brief determines whether interpolation from source to target is possible. |
467 |
Must be implemented by the actual Domain adapter |
468 |
*/ |
469 |
FINLEY_DLL_API |
470 |
virtual bool probeInterpolationACross(int functionSpaceType_source,const escript::AbstractDomain& targetDomain, int functionSpaceType_target) const; |
471 |
|
472 |
/** |
473 |
\brief |
474 |
copies the surface normals at data points into out. The actual function space to be considered |
475 |
is defined by out. out has to be defined on this. |
476 |
*/ |
477 |
FINLEY_DLL_API |
478 |
virtual void setToNormal(escript::Data& out) const; |
479 |
|
480 |
/** |
481 |
\brief |
482 |
copies the size of samples into out. The actual function space to be considered |
483 |
is defined by out. out has to be defined on this. |
484 |
*/ |
485 |
FINLEY_DLL_API |
486 |
virtual void setToSize(escript::Data& out) const; |
487 |
|
488 |
/** |
489 |
\brief |
490 |
copies the gradient of arg into grad. The actual function space to be considered |
491 |
for the gradient is defined by grad. arg and grad have to be defined on this. |
492 |
*/ |
493 |
FINLEY_DLL_API |
494 |
virtual void setToGradient(escript::Data& grad,const escript::Data& arg) const; |
495 |
|
496 |
/** |
497 |
\brief |
498 |
copies the integrals of the function defined by arg into integrals. |
499 |
arg has to be defined on this. |
500 |
*/ |
501 |
FINLEY_DLL_API |
502 |
virtual void setToIntegrals(std::vector<double>& integrals,const escript::Data& arg) const; |
503 |
|
504 |
/** |
505 |
\brief |
506 |
return the identifier of the matrix type to be used for the global stiffness matrix when a particular solver, package, perconditioner, |
507 |
and symmetric matrix is used. |
508 |
\param solver |
509 |
\param preconditioner |
510 |
\param package |
511 |
\param symmetry |
512 |
*/ |
513 |
FINLEY_DLL_API |
514 |
virtual int getSystemMatrixTypeId(const int solver, const int preconditioner, const int package, const bool symmetry) const; |
515 |
|
516 |
/** |
517 |
\brief |
518 |
return the identifier of the transport problem type to be used when a particular solver, perconditioner, package |
519 |
and symmetric matrix is used. |
520 |
\param solver |
521 |
\param preconditioner |
522 |
\param package |
523 |
\param symmetry |
524 |
*/ |
525 |
FINLEY_DLL_API |
526 |
virtual int getTransportTypeId(const int solver, const int preconditioner, const int package, const bool symmetry) const; |
527 |
|
528 |
/** |
529 |
\brief |
530 |
returns true if data on this domain and a function space of type functionSpaceCode has to |
531 |
considered as cell centered data. |
532 |
*/ |
533 |
FINLEY_DLL_API |
534 |
virtual bool isCellOriented(int functionSpaceCode) const; |
535 |
|
536 |
/** |
537 |
\brief |
538 |
Saves a dictonary of Data objects to an OpenDX input file. The keywords are used as identifier |
539 |
|
540 |
This has to be implemented by the actual Domain adapter. |
541 |
*/ |
542 |
FINLEY_DLL_API |
543 |
virtual void saveDX(const std::string& filename,const boost::python::dict& arg) const; |
544 |
|
545 |
|
546 |
/** |
547 |
\brief |
548 |
Saves a dictonary of Data objects to an VTK XML input file. The keywords are used as identifier |
549 |
|
550 |
This has to be implemented by the actual Domain adapter. |
551 |
*/ |
552 |
FINLEY_DLL_API |
553 |
virtual void saveVTK(const std::string& filename,const boost::python::dict& arg, const std::string& metadata, const std::string& metadata_schema) const; |
554 |
|
555 |
/** |
556 |
\brief |
557 |
returns the function space representation of the type functionSpaceCode on this domain |
558 |
as a vtkObject. |
559 |
*/ |
560 |
// vtkObject createVtkObject(int functionSpaceCode) const; |
561 |
|
562 |
/** |
563 |
\brief |
564 |
adds a PDE onto the stiffness matrix mat and a rhs |
565 |
*/ |
566 |
FINLEY_DLL_API |
567 |
virtual void addPDEToSystem( |
568 |
SystemMatrixAdapter& mat, escript::Data& rhs, |
569 |
const escript::Data& A, const escript::Data& B, const escript::Data& C, |
570 |
const escript::Data& D, const escript::Data& X, const escript::Data& Y, |
571 |
const escript::Data& d, const escript::Data& y, |
572 |
const escript::Data& d_contact, const escript::Data& y_contact) const; |
573 |
/** |
574 |
\brief |
575 |
adds a PDE onto the lumped stiffness matrix matrix |
576 |
*/ |
577 |
FINLEY_DLL_API |
578 |
virtual void addPDEToLumpedSystem( |
579 |
escript::Data& mat, |
580 |
const escript::Data& D, |
581 |
const escript::Data& d) const; |
582 |
|
583 |
/** |
584 |
\brief |
585 |
adds a PDE onto the stiffness matrix mat and a rhs |
586 |
*/ |
587 |
FINLEY_DLL_API |
588 |
virtual void addPDEToRHS(escript::Data& rhs, |
589 |
const escript::Data& X, const escript::Data& Y, |
590 |
const escript::Data& y, const escript::Data& y_contact) const; |
591 |
/** |
592 |
\brief |
593 |
adds a PDE onto a transport problem |
594 |
*/ |
595 |
|
596 |
FINLEY_DLL_API |
597 |
virtual void addPDEToTransportProblem( |
598 |
TransportProblemAdapter& tp, escript::Data& source, |
599 |
const escript::Data& M, |
600 |
const escript::Data& A, const escript::Data& B, const escript::Data& C,const escript::Data& D, |
601 |
const escript::Data& X,const escript::Data& Y, |
602 |
const escript::Data& d, const escript::Data& y, |
603 |
const escript::Data& d_contact,const escript::Data& y_contact) const; |
604 |
|
605 |
|
606 |
/** |
607 |
\brief |
608 |
creates a SystemMatrixAdapter stiffness matrix and initializes it with zeros: |
609 |
*/ |
610 |
FINLEY_DLL_API |
611 |
SystemMatrixAdapter newSystemMatrix( |
612 |
const int row_blocksize, |
613 |
const escript::FunctionSpace& row_functionspace, |
614 |
const int column_blocksize, |
615 |
const escript::FunctionSpace& column_functionspace, |
616 |
const int type) const; |
617 |
/** |
618 |
\brief |
619 |
creates a TransportProblemAdapter |
620 |
|
621 |
*/ |
622 |
|
623 |
FINLEY_DLL_API |
624 |
TransportProblemAdapter newTransportProblem( |
625 |
const double theta, |
626 |
const int blocksize, |
627 |
const escript::FunctionSpace& functionspace, |
628 |
const int type) const; |
629 |
|
630 |
/** |
631 |
\brief returns locations in the FEM nodes |
632 |
*/ |
633 |
FINLEY_DLL_API |
634 |
virtual escript::Data getX() const; |
635 |
|
636 |
/** |
637 |
\brief return boundary normals at the quadrature point on the face elements |
638 |
*/ |
639 |
FINLEY_DLL_API |
640 |
virtual escript::Data getNormal() const; |
641 |
|
642 |
/** |
643 |
\brief returns the element size |
644 |
*/ |
645 |
FINLEY_DLL_API |
646 |
virtual escript::Data getSize() const; |
647 |
|
648 |
/** |
649 |
\brief comparison operators |
650 |
*/ |
651 |
FINLEY_DLL_API |
652 |
virtual bool operator==(const escript::AbstractDomain& other) const; |
653 |
FINLEY_DLL_API |
654 |
virtual bool operator!=(const escript::AbstractDomain& other) const; |
655 |
|
656 |
/** |
657 |
\brief assigns new tag newTag to all samples of functionspace with a positive |
658 |
value of mask for any its sample point. |
659 |
|
660 |
*/ |
661 |
FINLEY_DLL_API |
662 |
virtual void setTags(const int functionSpaceType, const int newTag, const escript::Data& mask) const; |
663 |
|
664 |
/** |
665 |
\brief |
666 |
return the number of tags in use and a pointer to an array with the number of tags in use |
667 |
*/ |
668 |
FINLEY_DLL_API |
669 |
virtual int getNumberOfTagsInUse(int functionSpaceCode) const; |
670 |
|
671 |
FINLEY_DLL_API |
672 |
virtual const int* borrowListOfTagsInUse(int functionSpaceCode) const; |
673 |
|
674 |
|
675 |
/** |
676 |
\brief Checks if this domain allows tags for the specified functionSpaceCode. |
677 |
*/ |
678 |
FINLEY_DLL_API |
679 |
virtual |
680 |
bool canTag(int functionSpaceCode) const; |
681 |
|
682 |
|
683 |
protected: |
684 |
|
685 |
private: |
686 |
void extractArgsFromDict(const boost::python::dict& arg, int& numData, |
687 |
char**& names, escriptDataC*& data, |
688 |
escriptDataC**& dataPtr) const; |
689 |
|
690 |
// |
691 |
// pointer to the externally created finley mesh |
692 |
boost::shared_ptr<Finley_Mesh> m_finleyMesh; |
693 |
|
694 |
static FunctionSpaceNamesMapType m_functionSpaceTypeNames; |
695 |
|
696 |
}; |
697 |
|
698 |
} // end of namespace |
699 |
|
700 |
#endif |