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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1877 - (hide annotations)
Tue Oct 14 02:58:39 2008 UTC (14 years, 5 months ago) by ksteube
File MIME type: text/plain
File size: 17469 byte(s)
convection.py checkpointing uses mkdir/rmdir, and under MPI there
was a race condition.

mkdir needs to be run on only one CPU and then a barrier to prevent
working processors from using the directory before it exists.

Added methods domain.MPIBarrier and domain.onMasterProcessor() to
implement this technique.

A more general solution might be possible in the future.

1 ksteube 1312
2     /*******************************************************
3 ksteube 1811 *
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 ksteube 1312
14 ksteube 1811
15 jgs 149 #if !defined finley_MeshAdapter_20040526_H
16 jgs 82 #define finley_MeshAdapter_20040526_H
17 woo409 757 #include "system_dep.h"
18 jgs 82
19     extern "C" {
20 robwdcock 682 #include "../Mesh.h"
21     #include "../Finley.h"
22     #include "../Assemble.h"
23     #include "paso/SystemMatrix.h"
24 gross 1366 #include "paso/SolverFCT.h"
25 ksteube 1800 #include "paso/Paso_MPI.h"
26 jgs 82 }
27 jgs 472
28     #include "FinleyError.h"
29     #include "FinleyAdapterException.h"
30 jgs 149
31 jgs 480 #include "SystemMatrixAdapter.h"
32 gross 1366 #include "TransportProblemAdapter.h"
33 robwdcock 682 #include "escript/AbstractContinuousDomain.h"
34     #include "escript/FunctionSpace.h"
35     #include "escript/FunctionSpaceFactory.h"
36 jgs 472
37 jgs 82 #include <boost/shared_ptr.hpp>
38 jgs 153 #include <boost/python/dict.hpp>
39 jgs 472 #include <boost/python/extract.hpp>
40    
41 jgs 82 #include <map>
42     #include <vector>
43     #include <string>
44 jgs 472 #include <sstream>
45 jgs 82
46 jgs 480 //
47     // forward declarations
48     class Data;
49    
50     //using namespace escript;
51    
52 jgs 82 namespace finley {
53    
54 jgs 123 struct null_deleter
55     {
56     void operator()(void const *ptr) const
57     {
58     }
59     };
60    
61    
62 jgs 82 /**
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 jgs 480 class MeshAdapter : public escript::AbstractContinuousDomain {
73 jgs 82
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 gross 1062 static const int ReducedNodes;
82 jgs 82 static const int Elements;
83 gross 1059 static const int ReducedElements;
84 jgs 82 static const int FaceElements;
85 gross 1059 static const int ReducedFaceElements;
86 jgs 82 static const int Points;
87     static const int ContactElementsZero;
88 gross 1059 static const int ReducedContactElementsZero;
89 jgs 82 static const int ContactElementsOne;
90 gross 1059 static const int ReducedContactElementsOne;
91 jgs 82
92     /**
93     \brief
94     Constructor for MeshAdapter
95    
96     Description:
97     Constructor for MeshAdapter. The pointer passed to MeshAdapter
98 ksteube 1312 is deleted using a call to Finley_Mesh_free in the
99 jgs 82 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 ksteube 1312 Finley_Mesh_free in the MeshAdapter
108 jgs 82 destructor.
109     */
110 woo409 757 FINLEY_DLL_API
111 jgs 82 MeshAdapter(Finley_Mesh* finleyMesh=0);
112 jgs 149
113 jgs 82 /**
114     \brief
115     Copy constructor.
116     */
117 woo409 757 FINLEY_DLL_API
118 jgs 82 MeshAdapter(const MeshAdapter& in);
119 jgs 149
120 jgs 82 /**
121     \brief
122     Destructor for MeshAdapter. As specified in the constructor
123 ksteube 1312 this calls Finley_Mesh_free for the pointer given to the
124 jgs 82 constructor.
125     */
126 woo409 757 FINLEY_DLL_API
127 jgs 82 ~MeshAdapter();
128 jgs 149
129 jgs 82 /**
130     \brief
131 ksteube 1312 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 ksteube 1877 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     /**
160     \brief
161 jgs 82 return this as an AbstractContinuousDomain.
162     */
163     inline const AbstractContinuousDomain& asAbstractContinuousDomain() const
164     {
165     return *(static_cast<const AbstractContinuousDomain*>(this));
166     }
167    
168     /**
169     \brief
170     return this as an AbstractDomain.
171     */
172     inline const AbstractDomain& asAbstractDomain() const
173     {
174     return *(static_cast<const AbstractDomain*>(this));
175     }
176 jgs 149
177 jgs 82 /**
178     \brief
179     Write the current mesh to a file with the given name.
180     \param fileName Input - The name of the file to write to.
181     */
182 woo409 757 FINLEY_DLL_API
183 jgs 82 void write(const std::string& fileName) const;
184 jgs 149
185 jgs 82 /**
186     \brief
187 ksteube 1326 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 ksteube 1339 void Print_Mesh_Info(const bool) const;
192 ksteube 1326
193     /**
194     \brief
195 ksteube 1312 dumps the mesh to a file with the given name.
196     \param fileName Input - The name of the file
197     */
198     FINLEY_DLL_API
199     void dump(const std::string& fileName) const;
200    
201     /**
202     \brief
203 jgs 149 return the pointer to the underlying finley mesh structure
204 jgs 82 */
205 woo409 757 FINLEY_DLL_API
206 jgs 82 Finley_Mesh* getFinley_Mesh() const;
207 jgs 149
208 jgs 110 /**
209     \brief
210     Return the tag key for the given sample number.
211     \param functionSpaceType Input - The function space type.
212     \param sampleNo Input - The sample number.
213     */
214 woo409 757 FINLEY_DLL_API
215 jgs 110 int getTagFromSampleNo(int functionSpaceType, int sampleNo) const;
216 jgs 149
217 jgs 110 /**
218     \brief
219     Return the reference number of the given sample number.
220     \param functionSpaceType Input - The function space type.
221     */
222 woo409 757 FINLEY_DLL_API
223 gross 964 int* borrowSampleReferenceIDs(int functionSpaceType) const;
224 jgs 110
225     /**
226     \brief
227 jgs 82 Returns true if the given integer is a valid function space type
228     for this domain.
229     */
230 woo409 757 FINLEY_DLL_API
231 jgs 82 virtual bool isValidFunctionSpaceType(int functionSpaceType) const;
232 jgs 149
233 jgs 82 /**
234     \brief
235     Return a description for this domain
236     */
237 woo409 757 FINLEY_DLL_API
238 jgs 82 virtual std::string getDescription() const;
239 jgs 149
240 jgs 82 /**
241     \brief
242     Return a description for the given function space type code
243     */
244 woo409 757 FINLEY_DLL_API
245 jgs 82 virtual std::string functionSpaceTypeAsString(int functionSpaceType) const;
246 jgs 110
247 jgs 82 /**
248     \brief
249     Build the table of function space type names
250     */
251 woo409 757 FINLEY_DLL_API
252 jgs 82 void setFunctionSpaceTypeNames();
253 jgs 149
254 jgs 82 /**
255     \brief
256     Return a continuous FunctionSpace code
257     */
258 woo409 757 FINLEY_DLL_API
259 jgs 82 virtual int getContinuousFunctionCode() const;
260 jgs 149
261 jgs 82 /**
262     \brief
263 gross 1062 Return a continuous on reduced order nodes FunctionSpace code
264     */
265     FINLEY_DLL_API
266     virtual int getReducedContinuousFunctionCode() const;
267    
268     /**
269     \brief
270 gross 1059 Return a function FunctionSpace code
271 jgs 82 */
272 woo409 757 FINLEY_DLL_API
273 jgs 82 virtual int getFunctionCode() const;
274 jgs 149
275 jgs 82 /**
276     \brief
277 gross 1059 Return a function with reduced integration order FunctionSpace code
278     */
279     FINLEY_DLL_API
280     virtual int getReducedFunctionCode() const;
281    
282     /**
283     \brief
284 jgs 82 Return a function on boundary FunctionSpace code
285     */
286 woo409 757 FINLEY_DLL_API
287 jgs 82 virtual int getFunctionOnBoundaryCode() const;
288 jgs 149
289 jgs 82 /**
290     \brief
291 gross 1059 Return a function on boundary with reduced integration order FunctionSpace code
292     */
293     FINLEY_DLL_API
294     virtual int getReducedFunctionOnBoundaryCode() const;
295    
296     /**
297     \brief
298 jgs 82 Return a FunctionOnContactZero code
299     */
300 woo409 757 FINLEY_DLL_API
301 jgs 82 virtual int getFunctionOnContactZeroCode() const;
302 jgs 149
303 jgs 82 /**
304     \brief
305 gross 1059 Return a FunctionOnContactZero code with reduced integration order
306     */
307     FINLEY_DLL_API
308     virtual int getReducedFunctionOnContactZeroCode() const;
309    
310     /**
311     \brief
312 jgs 82 Return a FunctionOnContactOne code
313     */
314 woo409 757 FINLEY_DLL_API
315 jgs 82 virtual int getFunctionOnContactOneCode() const;
316 jgs 149
317 jgs 82 /**
318     \brief
319 gross 1059 Return a FunctionOnContactOne code with reduced integration order
320     */
321     FINLEY_DLL_API
322     virtual int getReducedFunctionOnContactOneCode() const;
323    
324     /**
325     \brief
326 jgs 82 Return a Solution code
327     */
328 woo409 757 FINLEY_DLL_API
329 jgs 82 virtual int getSolutionCode() const;
330 jgs 149
331 jgs 82 /**
332     \brief
333     Return a ReducedSolution code
334     */
335 woo409 757 FINLEY_DLL_API
336 jgs 82 virtual int getReducedSolutionCode() const;
337 jgs 149
338 jgs 82 /**
339     \brief
340     Return a DiracDeltaFunction code
341     */
342 woo409 757 FINLEY_DLL_API
343 jgs 82 virtual int getDiracDeltaFunctionCode() const;
344 jgs 149
345     /**
346 bcumming 782 5B
347 jgs 149 \brief
348     */
349 jgs 82 typedef std::map<int, std::string> FunctionSpaceNamesMapType;
350 jgs 149
351 jgs 82 /**
352     \brief
353     */
354 woo409 757 FINLEY_DLL_API
355 jgs 82 virtual int getDim() const;
356 jgs 149
357     /**
358 jgs 82 \brief
359 ksteube 1754 Return the number of data points summed across all MPI processes
360     */
361     FINLEY_DLL_API
362     virtual int getNumDataPointsGlobal() const;
363    
364     /**
365     \brief
366 jgs 82 Return the number of data points per sample, and the number of samples as a pair.
367 jgs 121 \param functionSpaceCode Input -
368 jgs 82 */
369 woo409 757 FINLEY_DLL_API
370 jgs 82 virtual std::pair<int,int> getDataShape(int functionSpaceCode) const;
371    
372     /**
373     \brief
374     copies the location of data points into arg. The domain of arg has to match this.
375     has to be implemented by the actual Domain adapter.
376     */
377 woo409 757 FINLEY_DLL_API
378 jgs 82 virtual void setToX(escript::Data& arg) const;
379 jgs 149
380 jgs 82 /**
381     \brief
382 gross 1044 sets a map from a clear tag name to a tag key
383     \param name Input - tag name.
384     \param tag Input - tag key.
385     */
386     FINLEY_DLL_API
387     virtual void setTagMap(const std::string& name, int tag);
388    
389     /**
390     \brief
391     Return the tag key for tag name.
392     \param name Input - tag name
393     */
394     FINLEY_DLL_API
395     virtual int getTag(const std::string& name) const;
396    
397     /**
398     \brief
399     Returns true if name is a defined tage name.
400     \param name Input - tag name to be checked.
401     */
402     FINLEY_DLL_API
403     virtual bool isValidTagName(const std::string& name) const;
404    
405     /**
406     \brief
407     Returns all tag names in a single string sperated by commas
408     */
409     FINLEY_DLL_API
410     virtual std::string showTagNames() const;
411    
412     /**
413     \brief
414 jgs 82 assigns new location to the domain
415     */
416 woo409 757 FINLEY_DLL_API
417 jgs 82 virtual void setNewX(const escript::Data& arg);
418 jgs 149
419 jgs 82 /**
420     \brief
421     interpolates data given on source onto target where source and target have to be given on the same domain.
422     */
423 woo409 757 FINLEY_DLL_API
424 jgs 82 virtual void interpolateOnDomain(escript::Data& target,const escript::Data& source) const;
425 woo409 757 FINLEY_DLL_API
426 jgs 82 virtual bool probeInterpolationOnDomain(int functionSpaceType_source,int functionSpaceType_target) const;
427 jgs 149
428 jgs 82 /**
429     \brief
430     interpolates data given on source onto target where source and target are given on different domains.
431     has to be implemented by the actual Domain adapter.
432     */
433 woo409 757 FINLEY_DLL_API
434 jgs 82 virtual void interpolateACross(escript::Data& target, const escript::Data& source) const;
435 woo409 757 FINLEY_DLL_API
436 jgs 82 virtual bool probeInterpolationACross(int functionSpaceType_source,const AbstractDomain& targetDomain, int functionSpaceType_target) const;
437 jgs 149
438 jgs 82 /**
439     \brief
440 jgs 149 copies the surface normals at data points into out. The actual function space to be considered
441 jgs 82 is defined by out. out has to be defined on this.
442     */
443 woo409 757 FINLEY_DLL_API
444 jgs 82 virtual void setToNormal(escript::Data& out) const;
445 jgs 149
446 jgs 82 /**
447     \brief
448     copies the size of samples into out. The actual function space to be considered
449     is defined by out. out has to be defined on this.
450     */
451 woo409 757 FINLEY_DLL_API
452 jgs 82 virtual void setToSize(escript::Data& out) const;
453    
454     /**
455     \brief
456     copies the gradient of arg into grad. The actual function space to be considered
457     for the gradient is defined by grad. arg and grad have to be defined on this.
458     */
459 woo409 757 FINLEY_DLL_API
460 jgs 82 virtual void setToGradient(escript::Data& grad,const escript::Data& arg) const;
461    
462     /**
463     \brief
464     copies the integrals of the function defined by arg into integrals.
465     arg has to be defined on this.
466     */
467 woo409 757 FINLEY_DLL_API
468 jgs 82 virtual void setToIntegrals(std::vector<double>& integrals,const escript::Data& arg) const;
469    
470 jgs 149 /**
471 jgs 102 \brief
472 gross 1859 return the identifier of the matrix type to be used for the global stiffness matrix when a particular solver, package, perconditioner,
473 jgs 102 and symmetric matrix is used.
474 gross 1859 \param precondioner
475 jgs 102 \param solver
476     \param symmetry
477     */
478 woo409 757 FINLEY_DLL_API
479 gross 1859 virtual int getSystemMatrixTypeId(const int solver, const int preconditioner, const int package, const bool symmetry) const;
480 jgs 102
481 jgs 82 /**
482     \brief
483 gross 1859 return the identifier of the transport problem type to be used when a particular solver, perconditioner, package
484     and symmetric matrix is used.
485     \param precondioner
486     \param solver
487     \param symmetry
488     */
489     FINLEY_DLL_API
490     virtual int getTransportTypeId(const int solver, const int preconditioner, const int package, const bool symmetry) const;
491    
492     /**
493     \brief
494 jgs 82 returns true if data on this domain and a function space of type functionSpaceCode has to
495     considered as cell centered data.
496     */
497 woo409 757 FINLEY_DLL_API
498 jgs 82 virtual bool isCellOriented(int functionSpaceCode) const;
499 jgs 149
500 jgs 82 /**
501     \brief
502 jgs 153 Saves a dictonary of Data objects to an OpenDX input file. The keywords are used as identifier
503    
504     This has to be implemented by the actual Domain adapter.
505 jgs 82 */
506 woo409 757 FINLEY_DLL_API
507 jgs 153 virtual void saveDX(const std::string& filename,const boost::python::dict& arg) const;
508 jgs 149
509 jgs 153
510 jgs 82 /**
511     \brief
512 jgs 153 Saves a dictonary of Data objects to an VTK XML input file. The keywords are used as identifier
513    
514     This has to be implemented by the actual Domain adapter.
515 jgs 110 */
516 woo409 757 FINLEY_DLL_API
517 jgs 153 virtual void saveVTK(const std::string& filename,const boost::python::dict& arg) const;
518 jgs 149
519 jgs 110 /**
520     \brief
521 jgs 82 returns the function space representation of the type functionSpaceCode on this domain
522     as a vtkObject.
523     */
524     // vtkObject createVtkObject(int functionSpaceCode) const;
525 jgs 149
526 jgs 82 /**
527     \brief
528     adds a PDE onto the stiffness matrix mat and a rhs
529     */
530 woo409 757 FINLEY_DLL_API
531 jgs 82 virtual void addPDEToSystem(
532     SystemMatrixAdapter& mat, escript::Data& rhs,
533     const escript::Data& A, const escript::Data& B, const escript::Data& C,
534 jgs 102 const escript::Data& D, const escript::Data& X, const escript::Data& Y,
535     const escript::Data& d, const escript::Data& y,
536     const escript::Data& d_contact, const escript::Data& y_contact) const;
537 gross 1204 /**
538     \brief
539     adds a PDE onto the lumped stiffness matrix matrix
540     */
541     FINLEY_DLL_API
542     virtual void addPDEToLumpedSystem(
543     escript::Data& mat,
544     const escript::Data& D,
545     const escript::Data& d) const;
546 jgs 149
547 jgs 82 /**
548     \brief
549 jgs 102 adds a PDE onto the stiffness matrix mat and a rhs
550 jgs 82 */
551 woo409 757 FINLEY_DLL_API
552 jgs 102 virtual void addPDEToRHS(escript::Data& rhs,
553     const escript::Data& X, const escript::Data& Y,
554     const escript::Data& y, const escript::Data& y_contact) const;
555 gross 1367 /**
556     \brief
557     adds a PDE onto a transport problem
558     */
559 jgs 149
560 gross 1367 FINLEY_DLL_API
561     virtual void addPDEToTransportProblem(
562     TransportProblemAdapter& tp, escript::Data& source,
563     const escript::Data& M,
564     const escript::Data& A, const escript::Data& B, const escript::Data& C,const escript::Data& D,
565     const escript::Data& X,const escript::Data& Y,
566     const escript::Data& d, const escript::Data& y,
567     const escript::Data& d_contact,const escript::Data& y_contact) const;
568    
569    
570 jgs 82 /**
571     \brief
572 gross 1366 creates a SystemMatrixAdapter stiffness matrix and initializes it with zeros:
573 jgs 82 */
574 woo409 757 FINLEY_DLL_API
575 jgs 82 SystemMatrixAdapter newSystemMatrix(
576     const int row_blocksize,
577     const escript::FunctionSpace& row_functionspace,
578     const int column_blocksize,
579     const escript::FunctionSpace& column_functionspace,
580 jgs 102 const int type) const;
581 gross 1366 /**
582     \brief
583     creates a TransportProblemAdapter
584 jgs 82
585 gross 1366 */
586    
587     FINLEY_DLL_API
588     TransportProblemAdapter newTransportProblem(
589     const double theta,
590     const int blocksize,
591     const escript::FunctionSpace& functionspace,
592     const int type) const;
593    
594 jgs 102 /**
595     \brief returns locations in the FEM nodes
596     */
597 woo409 757 FINLEY_DLL_API
598 jgs 102 virtual escript::Data getX() const;
599 jgs 149
600 jgs 102 /**
601     \brief return boundary normals at the quadrature point on the face elements
602     */
603 woo409 757 FINLEY_DLL_API
604 jgs 102 virtual escript::Data getNormal() const;
605 jgs 149
606 jgs 102 /**
607     \brief returns the element size
608     */
609 woo409 757 FINLEY_DLL_API
610 jgs 102 virtual escript::Data getSize() const;
611    
612 jgs 149 /**
613     \brief comparison operators
614     */
615 woo409 757 FINLEY_DLL_API
616 jgs 121 virtual bool operator==(const AbstractDomain& other) const;
617 woo409 757 FINLEY_DLL_API
618 jgs 121 virtual bool operator!=(const AbstractDomain& other) const;
619 jgs 102
620 gross 767 /**
621     \brief assigns new tag newTag to all samples of functionspace with a positive
622     value of mask for any its sample point.
623    
624     */
625     FINLEY_DLL_API
626     virtual void setTags(const int functionSpaceType, const int newTag, const escript::Data& mask) const;
627    
628 gross 1716 /**
629     \brief
630     return the number of tags in use and a pointer to an array with the number of tags in use
631     */
632     FINLEY_DLL_API
633     virtual int getNumberOfTagsInUse(int functionSpaceCode) const;
634    
635     FINLEY_DLL_API
636     virtual int* borrowListOfTagsInUse(int functionSpaceCode) const;
637    
638 jfenwick 1802
639     /**
640     \brief Checks if this domain allows tags for the specified functionSpaceCode.
641     */
642     FINLEY_DLL_API
643     virtual
644     bool canTag(int functionSpaceCode) const;
645    
646    
647 jgs 82 protected:
648    
649     private:
650 jgs 149
651 jgs 82 //
652     // pointer to the externally created finley mesh
653     boost::shared_ptr<Finley_Mesh> m_finleyMesh;
654    
655     static FunctionSpaceNamesMapType m_functionSpaceTypeNames;
656    
657     };
658    
659     } // end of namespace
660 jgs 149
661 jgs 82 #endif

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26