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

Contents of /trunk/finley/src/CPPAdapter/MeshAdapter.h

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1859 - (show annotations)
Wed Oct 8 03:03:37 2008 UTC (14 years, 5 months ago) by gross
File MIME type: text/plain
File size: 17187 byte(s)
first version of testing for transport solver.
1
2 /*******************************************************
3 *
4 * Copyright (c) 2003-2008 by University of Queensland
5 * Earth Systems Science Computational Center (ESSCC)
6 * http://www.uq.edu.au/esscc
7 *
8 * Primary Business: Queensland, Australia
9 * Licensed under the Open Software License version 3.0
10 * http://www.opensource.org/licenses/osl-3.0.php
11 *
12 *******************************************************/
13
14
15 #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 return this as an AbstractContinuousDomain.
146 */
147 inline const AbstractContinuousDomain& asAbstractContinuousDomain() const
148 {
149 return *(static_cast<const AbstractContinuousDomain*>(this));
150 }
151
152 /**
153 \brief
154 return this as an AbstractDomain.
155 */
156 inline const AbstractDomain& asAbstractDomain() const
157 {
158 return *(static_cast<const AbstractDomain*>(this));
159 }
160
161 /**
162 \brief
163 Write the current mesh to a file with the given name.
164 \param fileName Input - The name of the file to write to.
165 */
166 FINLEY_DLL_API
167 void write(const std::string& fileName) const;
168
169 /**
170 \brief
171 Write the current mesh to a file with the given name.
172 \param fileName Input - The name of the file to write to.
173 */
174 FINLEY_DLL_API
175 void Print_Mesh_Info(const bool) const;
176
177 /**
178 \brief
179 dumps the mesh to a file with the given name.
180 \param fileName Input - The name of the file
181 */
182 FINLEY_DLL_API
183 void dump(const std::string& fileName) const;
184
185 /**
186 \brief
187 return the pointer to the underlying finley mesh structure
188 */
189 FINLEY_DLL_API
190 Finley_Mesh* getFinley_Mesh() const;
191
192 /**
193 \brief
194 Return the tag key for the given sample number.
195 \param functionSpaceType Input - The function space type.
196 \param sampleNo Input - The sample number.
197 */
198 FINLEY_DLL_API
199 int getTagFromSampleNo(int functionSpaceType, int sampleNo) const;
200
201 /**
202 \brief
203 Return the reference number of the given sample number.
204 \param functionSpaceType Input - The function space type.
205 */
206 FINLEY_DLL_API
207 int* borrowSampleReferenceIDs(int functionSpaceType) const;
208
209 /**
210 \brief
211 Returns true if the given integer is a valid function space type
212 for this domain.
213 */
214 FINLEY_DLL_API
215 virtual bool isValidFunctionSpaceType(int functionSpaceType) const;
216
217 /**
218 \brief
219 Return a description for this domain
220 */
221 FINLEY_DLL_API
222 virtual std::string getDescription() const;
223
224 /**
225 \brief
226 Return a description for the given function space type code
227 */
228 FINLEY_DLL_API
229 virtual std::string functionSpaceTypeAsString(int functionSpaceType) const;
230
231 /**
232 \brief
233 Build the table of function space type names
234 */
235 FINLEY_DLL_API
236 void setFunctionSpaceTypeNames();
237
238 /**
239 \brief
240 Return a continuous FunctionSpace code
241 */
242 FINLEY_DLL_API
243 virtual int getContinuousFunctionCode() const;
244
245 /**
246 \brief
247 Return a continuous on reduced order nodes FunctionSpace code
248 */
249 FINLEY_DLL_API
250 virtual int getReducedContinuousFunctionCode() const;
251
252 /**
253 \brief
254 Return a function FunctionSpace code
255 */
256 FINLEY_DLL_API
257 virtual int getFunctionCode() const;
258
259 /**
260 \brief
261 Return a function with reduced integration order FunctionSpace code
262 */
263 FINLEY_DLL_API
264 virtual int getReducedFunctionCode() const;
265
266 /**
267 \brief
268 Return a function on boundary FunctionSpace code
269 */
270 FINLEY_DLL_API
271 virtual int getFunctionOnBoundaryCode() const;
272
273 /**
274 \brief
275 Return a function on boundary with reduced integration order FunctionSpace code
276 */
277 FINLEY_DLL_API
278 virtual int getReducedFunctionOnBoundaryCode() const;
279
280 /**
281 \brief
282 Return a FunctionOnContactZero code
283 */
284 FINLEY_DLL_API
285 virtual int getFunctionOnContactZeroCode() const;
286
287 /**
288 \brief
289 Return a FunctionOnContactZero code with reduced integration order
290 */
291 FINLEY_DLL_API
292 virtual int getReducedFunctionOnContactZeroCode() const;
293
294 /**
295 \brief
296 Return a FunctionOnContactOne code
297 */
298 FINLEY_DLL_API
299 virtual int getFunctionOnContactOneCode() const;
300
301 /**
302 \brief
303 Return a FunctionOnContactOne code with reduced integration order
304 */
305 FINLEY_DLL_API
306 virtual int getReducedFunctionOnContactOneCode() const;
307
308 /**
309 \brief
310 Return a Solution code
311 */
312 FINLEY_DLL_API
313 virtual int getSolutionCode() const;
314
315 /**
316 \brief
317 Return a ReducedSolution code
318 */
319 FINLEY_DLL_API
320 virtual int getReducedSolutionCode() const;
321
322 /**
323 \brief
324 Return a DiracDeltaFunction code
325 */
326 FINLEY_DLL_API
327 virtual int getDiracDeltaFunctionCode() const;
328
329 /**
330 5B
331 \brief
332 */
333 typedef std::map<int, std::string> FunctionSpaceNamesMapType;
334
335 /**
336 \brief
337 */
338 FINLEY_DLL_API
339 virtual int getDim() const;
340
341 /**
342 \brief
343 Return the number of data points summed across all MPI processes
344 */
345 FINLEY_DLL_API
346 virtual int getNumDataPointsGlobal() const;
347
348 /**
349 \brief
350 Return the number of data points per sample, and the number of samples as a pair.
351 \param functionSpaceCode Input -
352 */
353 FINLEY_DLL_API
354 virtual std::pair<int,int> getDataShape(int functionSpaceCode) const;
355
356 /**
357 \brief
358 copies the location of data points into arg. The domain of arg has to match this.
359 has to be implemented by the actual Domain adapter.
360 */
361 FINLEY_DLL_API
362 virtual void setToX(escript::Data& arg) const;
363
364 /**
365 \brief
366 sets a map from a clear tag name to a tag key
367 \param name Input - tag name.
368 \param tag Input - tag key.
369 */
370 FINLEY_DLL_API
371 virtual void setTagMap(const std::string& name, int tag);
372
373 /**
374 \brief
375 Return the tag key for tag name.
376 \param name Input - tag name
377 */
378 FINLEY_DLL_API
379 virtual int getTag(const std::string& name) const;
380
381 /**
382 \brief
383 Returns true if name is a defined tage name.
384 \param name Input - tag name to be checked.
385 */
386 FINLEY_DLL_API
387 virtual bool isValidTagName(const std::string& name) const;
388
389 /**
390 \brief
391 Returns all tag names in a single string sperated by commas
392 */
393 FINLEY_DLL_API
394 virtual std::string showTagNames() const;
395
396 /**
397 \brief
398 assigns new location to the domain
399 */
400 FINLEY_DLL_API
401 virtual void setNewX(const escript::Data& arg);
402
403 /**
404 \brief
405 interpolates data given on source onto target where source and target have to be given on the same domain.
406 */
407 FINLEY_DLL_API
408 virtual void interpolateOnDomain(escript::Data& target,const escript::Data& source) const;
409 FINLEY_DLL_API
410 virtual bool probeInterpolationOnDomain(int functionSpaceType_source,int functionSpaceType_target) const;
411
412 /**
413 \brief
414 interpolates data given on source onto target where source and target are given on different domains.
415 has to be implemented by the actual Domain adapter.
416 */
417 FINLEY_DLL_API
418 virtual void interpolateACross(escript::Data& target, const escript::Data& source) const;
419 FINLEY_DLL_API
420 virtual bool probeInterpolationACross(int functionSpaceType_source,const AbstractDomain& targetDomain, int functionSpaceType_target) const;
421
422 /**
423 \brief
424 copies the surface normals at data points into out. The actual function space to be considered
425 is defined by out. out has to be defined on this.
426 */
427 FINLEY_DLL_API
428 virtual void setToNormal(escript::Data& out) const;
429
430 /**
431 \brief
432 copies the size of samples into out. The actual function space to be considered
433 is defined by out. out has to be defined on this.
434 */
435 FINLEY_DLL_API
436 virtual void setToSize(escript::Data& out) const;
437
438 /**
439 \brief
440 copies the gradient of arg into grad. The actual function space to be considered
441 for the gradient is defined by grad. arg and grad have to be defined on this.
442 */
443 FINLEY_DLL_API
444 virtual void setToGradient(escript::Data& grad,const escript::Data& arg) const;
445
446 /**
447 \brief
448 copies the integrals of the function defined by arg into integrals.
449 arg has to be defined on this.
450 */
451 FINLEY_DLL_API
452 virtual void setToIntegrals(std::vector<double>& integrals,const escript::Data& arg) const;
453
454 /**
455 \brief
456 return the identifier of the matrix type to be used for the global stiffness matrix when a particular solver, package, perconditioner,
457 and symmetric matrix is used.
458 \param precondioner
459 \param solver
460 \param symmetry
461 */
462 FINLEY_DLL_API
463 virtual int getSystemMatrixTypeId(const int solver, const int preconditioner, const int package, const bool symmetry) const;
464
465 /**
466 \brief
467 return the identifier of the transport problem type to be used when a particular solver, perconditioner, package
468 and symmetric matrix is used.
469 \param precondioner
470 \param solver
471 \param symmetry
472 */
473 FINLEY_DLL_API
474 virtual int getTransportTypeId(const int solver, const int preconditioner, const int package, const bool symmetry) const;
475
476 /**
477 \brief
478 returns true if data on this domain and a function space of type functionSpaceCode has to
479 considered as cell centered data.
480 */
481 FINLEY_DLL_API
482 virtual bool isCellOriented(int functionSpaceCode) const;
483
484 /**
485 \brief
486 Saves a dictonary of Data objects to an OpenDX input file. The keywords are used as identifier
487
488 This has to be implemented by the actual Domain adapter.
489 */
490 FINLEY_DLL_API
491 virtual void saveDX(const std::string& filename,const boost::python::dict& arg) const;
492
493
494 /**
495 \brief
496 Saves a dictonary of Data objects to an VTK XML input file. The keywords are used as identifier
497
498 This has to be implemented by the actual Domain adapter.
499 */
500 FINLEY_DLL_API
501 virtual void saveVTK(const std::string& filename,const boost::python::dict& arg) const;
502
503 /**
504 \brief
505 returns the function space representation of the type functionSpaceCode on this domain
506 as a vtkObject.
507 */
508 // vtkObject createVtkObject(int functionSpaceCode) const;
509
510 /**
511 \brief
512 adds a PDE onto the stiffness matrix mat and a rhs
513 */
514 FINLEY_DLL_API
515 virtual void addPDEToSystem(
516 SystemMatrixAdapter& mat, escript::Data& rhs,
517 const escript::Data& A, const escript::Data& B, const escript::Data& C,
518 const escript::Data& D, const escript::Data& X, const escript::Data& Y,
519 const escript::Data& d, const escript::Data& y,
520 const escript::Data& d_contact, const escript::Data& y_contact) const;
521 /**
522 \brief
523 adds a PDE onto the lumped stiffness matrix matrix
524 */
525 FINLEY_DLL_API
526 virtual void addPDEToLumpedSystem(
527 escript::Data& mat,
528 const escript::Data& D,
529 const escript::Data& d) const;
530
531 /**
532 \brief
533 adds a PDE onto the stiffness matrix mat and a rhs
534 */
535 FINLEY_DLL_API
536 virtual void addPDEToRHS(escript::Data& rhs,
537 const escript::Data& X, const escript::Data& Y,
538 const escript::Data& y, const escript::Data& y_contact) const;
539 /**
540 \brief
541 adds a PDE onto a transport problem
542 */
543
544 FINLEY_DLL_API
545 virtual void addPDEToTransportProblem(
546 TransportProblemAdapter& tp, escript::Data& source,
547 const escript::Data& M,
548 const escript::Data& A, const escript::Data& B, const escript::Data& C,const escript::Data& D,
549 const escript::Data& X,const escript::Data& Y,
550 const escript::Data& d, const escript::Data& y,
551 const escript::Data& d_contact,const escript::Data& y_contact) const;
552
553
554 /**
555 \brief
556 creates a SystemMatrixAdapter stiffness matrix and initializes it with zeros:
557 */
558 FINLEY_DLL_API
559 SystemMatrixAdapter newSystemMatrix(
560 const int row_blocksize,
561 const escript::FunctionSpace& row_functionspace,
562 const int column_blocksize,
563 const escript::FunctionSpace& column_functionspace,
564 const int type) const;
565 /**
566 \brief
567 creates a TransportProblemAdapter
568
569 */
570
571 FINLEY_DLL_API
572 TransportProblemAdapter newTransportProblem(
573 const double theta,
574 const int blocksize,
575 const escript::FunctionSpace& functionspace,
576 const int type) const;
577
578 /**
579 \brief returns locations in the FEM nodes
580 */
581 FINLEY_DLL_API
582 virtual escript::Data getX() const;
583
584 /**
585 \brief return boundary normals at the quadrature point on the face elements
586 */
587 FINLEY_DLL_API
588 virtual escript::Data getNormal() const;
589
590 /**
591 \brief returns the element size
592 */
593 FINLEY_DLL_API
594 virtual escript::Data getSize() const;
595
596 /**
597 \brief comparison operators
598 */
599 FINLEY_DLL_API
600 virtual bool operator==(const AbstractDomain& other) const;
601 FINLEY_DLL_API
602 virtual bool operator!=(const AbstractDomain& other) const;
603
604 /**
605 \brief assigns new tag newTag to all samples of functionspace with a positive
606 value of mask for any its sample point.
607
608 */
609 FINLEY_DLL_API
610 virtual void setTags(const int functionSpaceType, const int newTag, const escript::Data& mask) const;
611
612 /**
613 \brief
614 return the number of tags in use and a pointer to an array with the number of tags in use
615 */
616 FINLEY_DLL_API
617 virtual int getNumberOfTagsInUse(int functionSpaceCode) const;
618
619 FINLEY_DLL_API
620 virtual int* borrowListOfTagsInUse(int functionSpaceCode) const;
621
622
623 /**
624 \brief Checks if this domain allows tags for the specified functionSpaceCode.
625 */
626 FINLEY_DLL_API
627 virtual
628 bool canTag(int functionSpaceCode) const;
629
630
631 protected:
632
633 private:
634
635 //
636 // pointer to the externally created finley mesh
637 boost::shared_ptr<Finley_Mesh> m_finleyMesh;
638
639 static FunctionSpaceNamesMapType m_functionSpaceTypeNames;
640
641 };
642
643 } // end of namespace
644
645 #endif

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26