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 |
#if !defined finley_MeshAdapter_20040526_H |
17 |
#define finley_MeshAdapter_20040526_H |
18 |
|
19 |
extern "C" { |
20 |
#include "Mesh.h" |
21 |
#include "Finley.h" |
22 |
#include "Assemble.h" |
23 |
#include "Finley.h" |
24 |
#include "SystemMatrix.h" |
25 |
} |
26 |
|
27 |
#include "FinleyError.h" |
28 |
#include "FinleyAdapterException.h" |
29 |
|
30 |
#include "SystemMatrixAdapter.h" |
31 |
#include "AbstractContinuousDomain.h" |
32 |
#include "FunctionSpace.h" |
33 |
#include "FunctionSpaceFactory.h" |
34 |
|
35 |
#include <boost/shared_ptr.hpp> |
36 |
#include <boost/python/dict.hpp> |
37 |
#include <boost/python/extract.hpp> |
38 |
|
39 |
#include <map> |
40 |
#include <vector> |
41 |
#include <string> |
42 |
#include <sstream> |
43 |
|
44 |
// |
45 |
// forward declarations |
46 |
class Data; |
47 |
|
48 |
//using namespace escript; |
49 |
|
50 |
namespace finley { |
51 |
|
52 |
struct null_deleter |
53 |
{ |
54 |
void operator()(void const *ptr) const |
55 |
{ |
56 |
} |
57 |
}; |
58 |
|
59 |
|
60 |
/** |
61 |
\brief |
62 |
MeshAdapter implements the AbstractContinuousDomain |
63 |
interface for the Finley library. |
64 |
|
65 |
Description: |
66 |
MeshAdapter implements the AbstractContinuousDomain |
67 |
interface for the Finley library. |
68 |
*/ |
69 |
|
70 |
class MeshAdapter : public escript::AbstractContinuousDomain { |
71 |
|
72 |
public: |
73 |
|
74 |
// |
75 |
// Codes for function space types supported |
76 |
static const int DegreesOfFreedom; |
77 |
static const int ReducedDegreesOfFreedom; |
78 |
static const int Nodes; |
79 |
static const int Elements; |
80 |
static const int FaceElements; |
81 |
static const int Points; |
82 |
static const int ContactElementsZero; |
83 |
static const int ContactElementsOne; |
84 |
|
85 |
/** |
86 |
\brief |
87 |
Constructor for MeshAdapter |
88 |
|
89 |
Description: |
90 |
Constructor for MeshAdapter. The pointer passed to MeshAdapter |
91 |
is deleted using a call to Finley_Mesh_deallocate in the |
92 |
MeshAdapter destructor. |
93 |
|
94 |
Throws: |
95 |
May throw an exception derived from EsysException |
96 |
|
97 |
\param finleyMesh Input - A pointer to the externally constructed |
98 |
finley mesh.The pointer passed to MeshAdapter |
99 |
is deleted using a call to |
100 |
Finley_Mesh_deallocate in the MeshAdapter |
101 |
destructor. |
102 |
*/ |
103 |
MeshAdapter(Finley_Mesh* finleyMesh=0); |
104 |
|
105 |
/** |
106 |
\brief |
107 |
Copy constructor. |
108 |
*/ |
109 |
MeshAdapter(const MeshAdapter& in); |
110 |
|
111 |
/** |
112 |
\brief |
113 |
Destructor for MeshAdapter. As specified in the constructor |
114 |
this calls Finley_Mesh_deallocate for the pointer given to the |
115 |
constructor. |
116 |
*/ |
117 |
~MeshAdapter(); |
118 |
|
119 |
/** |
120 |
\brief |
121 |
return this as an AbstractContinuousDomain. |
122 |
*/ |
123 |
inline const AbstractContinuousDomain& asAbstractContinuousDomain() const |
124 |
{ |
125 |
return *(static_cast<const AbstractContinuousDomain*>(this)); |
126 |
} |
127 |
|
128 |
/** |
129 |
\brief |
130 |
return this as an AbstractDomain. |
131 |
*/ |
132 |
inline const AbstractDomain& asAbstractDomain() const |
133 |
{ |
134 |
return *(static_cast<const AbstractDomain*>(this)); |
135 |
} |
136 |
|
137 |
/** |
138 |
\brief |
139 |
Write the current mesh to a file with the given name. |
140 |
\param fileName Input - The name of the file to write to. |
141 |
*/ |
142 |
void write(const std::string& fileName) const; |
143 |
|
144 |
/** |
145 |
\brief |
146 |
return the pointer to the underlying finley mesh structure |
147 |
*/ |
148 |
Finley_Mesh* getFinley_Mesh() const; |
149 |
|
150 |
/** |
151 |
\brief |
152 |
Return the tag key for the given sample number. |
153 |
\param functionSpaceType Input - The function space type. |
154 |
\param sampleNo Input - The sample number. |
155 |
*/ |
156 |
int getTagFromSampleNo(int functionSpaceType, int sampleNo) const; |
157 |
|
158 |
/** |
159 |
\brief |
160 |
Return the reference number of the given sample number. |
161 |
\param functionSpaceType Input - The function space type. |
162 |
\param sampleNo Input - The sample number. |
163 |
*/ |
164 |
int getReferenceNoFromSampleNo(int functionSpaceType, int sampleNo) const; |
165 |
|
166 |
/** |
167 |
\brief |
168 |
Returns true if the given integer is a valid function space type |
169 |
for this domain. |
170 |
*/ |
171 |
virtual bool isValidFunctionSpaceType(int functionSpaceType) const; |
172 |
|
173 |
/** |
174 |
\brief |
175 |
Return a description for this domain |
176 |
*/ |
177 |
virtual std::string getDescription() const; |
178 |
|
179 |
/** |
180 |
\brief |
181 |
Return a description for the given function space type code |
182 |
*/ |
183 |
virtual std::string functionSpaceTypeAsString(int functionSpaceType) const; |
184 |
|
185 |
/** |
186 |
\brief |
187 |
Build the table of function space type names |
188 |
*/ |
189 |
void setFunctionSpaceTypeNames(); |
190 |
|
191 |
/** |
192 |
\brief |
193 |
Return a continuous FunctionSpace code |
194 |
*/ |
195 |
virtual int getContinuousFunctionCode() const; |
196 |
|
197 |
/** |
198 |
\brief |
199 |
Return a functon FunctionSpace code |
200 |
*/ |
201 |
virtual int getFunctionCode() const; |
202 |
|
203 |
/** |
204 |
\brief |
205 |
Return a function on boundary FunctionSpace code |
206 |
*/ |
207 |
virtual int getFunctionOnBoundaryCode() const; |
208 |
|
209 |
/** |
210 |
\brief |
211 |
Return a FunctionOnContactZero code |
212 |
*/ |
213 |
virtual int getFunctionOnContactZeroCode() const; |
214 |
|
215 |
/** |
216 |
\brief |
217 |
Return a FunctionOnContactOne code |
218 |
*/ |
219 |
virtual int getFunctionOnContactOneCode() const; |
220 |
|
221 |
/** |
222 |
\brief |
223 |
Return a Solution code |
224 |
*/ |
225 |
virtual int getSolutionCode() const; |
226 |
|
227 |
/** |
228 |
\brief |
229 |
Return a ReducedSolution code |
230 |
*/ |
231 |
virtual int getReducedSolutionCode() const; |
232 |
|
233 |
/** |
234 |
\brief |
235 |
Return a DiracDeltaFunction code |
236 |
*/ |
237 |
virtual int getDiracDeltaFunctionCode() const; |
238 |
|
239 |
/** |
240 |
\brief |
241 |
*/ |
242 |
typedef std::map<int, std::string> FunctionSpaceNamesMapType; |
243 |
|
244 |
/** |
245 |
\brief |
246 |
*/ |
247 |
virtual int getDim() const; |
248 |
|
249 |
/** |
250 |
\brief |
251 |
Return the number of data points per sample, and the number of samples as a pair. |
252 |
\param functionSpaceCode Input - |
253 |
*/ |
254 |
virtual std::pair<int,int> getDataShape(int functionSpaceCode) const; |
255 |
|
256 |
/** |
257 |
\brief |
258 |
copies the location of data points into arg. The domain of arg has to match this. |
259 |
has to be implemented by the actual Domain adapter. |
260 |
*/ |
261 |
virtual void setToX(escript::Data& arg) const; |
262 |
|
263 |
/** |
264 |
\brief |
265 |
assigns new location to the domain |
266 |
*/ |
267 |
virtual void setNewX(const escript::Data& arg); |
268 |
|
269 |
/** |
270 |
\brief |
271 |
interpolates data given on source onto target where source and target have to be given on the same domain. |
272 |
*/ |
273 |
virtual void interpolateOnDomain(escript::Data& target,const escript::Data& source) const; |
274 |
virtual bool probeInterpolationOnDomain(int functionSpaceType_source,int functionSpaceType_target) const; |
275 |
|
276 |
/** |
277 |
\brief |
278 |
interpolates data given on source onto target where source and target are given on different domains. |
279 |
has to be implemented by the actual Domain adapter. |
280 |
*/ |
281 |
virtual void interpolateACross(escript::Data& target, const escript::Data& source) const; |
282 |
virtual bool probeInterpolationACross(int functionSpaceType_source,const AbstractDomain& targetDomain, int functionSpaceType_target) const; |
283 |
|
284 |
/** |
285 |
\brief |
286 |
copies the surface normals at data points into out. The actual function space to be considered |
287 |
is defined by out. out has to be defined on this. |
288 |
*/ |
289 |
virtual void setToNormal(escript::Data& out) const; |
290 |
|
291 |
/** |
292 |
\brief |
293 |
copies the size of samples into out. The actual function space to be considered |
294 |
is defined by out. out has to be defined on this. |
295 |
*/ |
296 |
virtual void setToSize(escript::Data& out) const; |
297 |
|
298 |
/** |
299 |
\brief |
300 |
copies the gradient of arg into grad. The actual function space to be considered |
301 |
for the gradient is defined by grad. arg and grad have to be defined on this. |
302 |
*/ |
303 |
virtual void setToGradient(escript::Data& grad,const escript::Data& arg) const; |
304 |
|
305 |
/** |
306 |
\brief |
307 |
copies the integrals of the function defined by arg into integrals. |
308 |
arg has to be defined on this. |
309 |
*/ |
310 |
virtual void setToIntegrals(std::vector<double>& integrals,const escript::Data& arg) const; |
311 |
|
312 |
/** |
313 |
\brief |
314 |
return the identifier of the matrix type to be used for the global stiffness matrix when a particular solver, package |
315 |
and symmetric matrix is used. |
316 |
\param solver |
317 |
\param symmetry |
318 |
*/ |
319 |
virtual int getSystemMatrixTypeId(const int solver, const int package, const bool symmetry) const; |
320 |
|
321 |
/** |
322 |
\brief |
323 |
returns true if data on this domain and a function space of type functionSpaceCode has to |
324 |
considered as cell centered data. |
325 |
*/ |
326 |
virtual bool isCellOriented(int functionSpaceCode) const; |
327 |
|
328 |
/** |
329 |
\brief |
330 |
Saves a dictonary of Data objects to an OpenDX input file. The keywords are used as identifier |
331 |
|
332 |
This has to be implemented by the actual Domain adapter. |
333 |
*/ |
334 |
virtual void saveDX(const std::string& filename,const boost::python::dict& arg) const; |
335 |
|
336 |
|
337 |
/** |
338 |
\brief |
339 |
Saves a dictonary of Data objects to an VTK XML input file. The keywords are used as identifier |
340 |
|
341 |
This has to be implemented by the actual Domain adapter. |
342 |
*/ |
343 |
virtual void saveVTK(const std::string& filename,const boost::python::dict& arg) const; |
344 |
|
345 |
/** |
346 |
\brief |
347 |
returns the function space representation of the type functionSpaceCode on this domain |
348 |
as a vtkObject. |
349 |
*/ |
350 |
// vtkObject createVtkObject(int functionSpaceCode) const; |
351 |
|
352 |
/** |
353 |
\brief |
354 |
adds a PDE onto the stiffness matrix mat and a rhs |
355 |
*/ |
356 |
virtual void addPDEToSystem( |
357 |
SystemMatrixAdapter& mat, escript::Data& rhs, |
358 |
const escript::Data& A, const escript::Data& B, const escript::Data& C, |
359 |
const escript::Data& D, const escript::Data& X, const escript::Data& Y, |
360 |
const escript::Data& d, const escript::Data& y, |
361 |
const escript::Data& d_contact, const escript::Data& y_contact) const; |
362 |
|
363 |
/** |
364 |
\brief |
365 |
adds a PDE onto the stiffness matrix mat and a rhs |
366 |
*/ |
367 |
virtual void addPDEToRHS(escript::Data& rhs, |
368 |
const escript::Data& X, const escript::Data& Y, |
369 |
const escript::Data& y, const escript::Data& y_contact) const; |
370 |
|
371 |
/** |
372 |
\brief |
373 |
creates a SystemMatrixAdapter stiffness matrix an initializes it with zeros: |
374 |
*/ |
375 |
SystemMatrixAdapter newSystemMatrix( |
376 |
const int row_blocksize, |
377 |
const escript::FunctionSpace& row_functionspace, |
378 |
const int column_blocksize, |
379 |
const escript::FunctionSpace& column_functionspace, |
380 |
const int type) const; |
381 |
|
382 |
/** |
383 |
\brief returns locations in the FEM nodes |
384 |
*/ |
385 |
virtual escript::Data getX() const; |
386 |
|
387 |
/** |
388 |
\brief return boundary normals at the quadrature point on the face elements |
389 |
*/ |
390 |
virtual escript::Data getNormal() const; |
391 |
|
392 |
/** |
393 |
\brief returns the element size |
394 |
*/ |
395 |
virtual escript::Data getSize() const; |
396 |
|
397 |
/** |
398 |
\brief comparison operators |
399 |
*/ |
400 |
virtual bool operator==(const AbstractDomain& other) const; |
401 |
virtual bool operator!=(const AbstractDomain& other) const; |
402 |
|
403 |
protected: |
404 |
|
405 |
private: |
406 |
|
407 |
// |
408 |
// pointer to the externally created finley mesh |
409 |
boost::shared_ptr<Finley_Mesh> m_finleyMesh; |
410 |
|
411 |
static FunctionSpaceNamesMapType m_functionSpaceTypeNames; |
412 |
|
413 |
}; |
414 |
|
415 |
} // end of namespace |
416 |
|
417 |
#endif |