AFEPack
MPI_Migration.h
浏览该文件的文档。
00001 
00011 #ifndef __MPI_Migration_h__
00012 #define __MPI_Migration_h__
00013 
00014 #include <mpi.h>
00015 #include <boost/program_options.hpp>
00016 #include <AFEPack/Migration.details.h>
00017 #include <AFEPack/mpi/MPI_HGeometry.h>
00018 #include <AFEPack/FEMSpace.h>
00019 
00020 namespace Migration {
00021 
00022   namespace details {
00024 
00027     template <int DOW> void
00028       clear_data_buffer(HGeometry<0,DOW>& geo) {
00029       geo.buffer.clear();
00030     }
00031     template <int DOW> void
00032       clear_data_buffer(HGeometry<1,DOW>& geo) {
00033       geo.buffer.clear();
00034       for (u_int i = 0;i < geo.n_vertex;++ i) {
00035         clear_data_buffer(*geo.vertex[i]);
00036       }
00037       if (geo.isRefined()) {
00038         for (u_int i = 0;i < geo.n_child;++ i) {
00039           clear_data_buffer(*geo.child[i]);
00040         }
00041       }
00042     }
00043     template <class GEO> void
00044       clear_data_buffer(GEO& geo) {
00045       geo.buffer.clear();
00046       for (u_int i = 0;i < geo.n_boundary;++ i) {
00047         clear_data_buffer(*geo.boundary[i]);
00048       }
00049       if (geo.isRefined()) {
00050         for (u_int i = 0;i < geo.n_child;++ i) {
00051           clear_data_buffer(*geo.child[i]);
00052         }
00053       }
00054     }
00056   }
00057 
00063   template <class HTREE> void
00064     clear_data_buffer(HTREE& tree) {
00065     typename HTREE::RootIterator 
00066       the_ele = tree.beginRootElement(),
00067       end_ele = tree.endRootElement();
00068     for (;the_ele != end_ele;++ the_ele) {
00069       details::clear_data_buffer(*the_ele);
00070     }
00071   }
00072 
00073   void initialize(MPI_Comm comm); 
00074   data_id_t register_data_name(const data_name_t& dn); 
00075   void load_config(const std::string& filename); 
00076   void save_config(const std::string& filename); 
00077 
00085   template <class FUNC>
00086     void export_fe_func(const FUNC& fun,
00087                         const data_id_t& data_id) {
00088     typedef FUNC fe_func_t;
00089     typedef typename fe_func_t::fe_space_t fe_space_t;
00090     typedef typename fe_space_t::dof_info_t dof_info_t;
00091     typedef RegularMesh<fe_space_t::dim,fe_space_t::dow> mesh_t;
00092     const fe_space_t& sp = fun.femSpace();
00093     const mesh_t& mesh = dynamic_cast<const mesh_t&>(sp.mesh());
00094     AFEPack::ostream<> os;
00095     u_int n_dof = sp.n_dof();
00096     for (u_int i = 0;i < n_dof;++ i) {
00097       const dof_info_t& di = sp.dofInfo(i);
00098       const DOFIndex& dof_idx = sp.dofIndex(i);
00099       get_export_stream(mesh,
00100                         data_id,
00101                         dof_idx.dimension,
00102                         dof_idx.geometry_index,
00103                         os);
00104       os << di.identity
00105          << di.interp_point 
00106          << fun(i);
00107     }
00108   }
00118   template <class FUNC>
00119     void export_fe_func(const std::vector<FUNC>& fun,
00120                         const data_id_t& data_id) {
00121     typedef FUNC fe_func_t;
00122     typedef typename fe_func_t::fe_space_t fe_space_t;
00123     typedef typename fe_space_t::dof_info_t dof_info_t;
00124     typedef RegularMesh<fe_space_t::dim,fe_space_t::dow> mesh_t;
00125     const fe_space_t& sp = fun[0].femSpace();
00126     const mesh_t& mesh = dynamic_cast<const mesh_t&>(sp.mesh());
00127     AFEPack::ostream<> os;
00128     u_int n_dof = sp.n_dof();
00129     for (u_int i = 0;i < n_dof;++ i) {
00130       const dof_info_t& di = sp.dofInfo(i);
00131       const DOFIndex& dof_idx = sp.dofIndex(i);
00132       get_export_stream(mesh,
00133                         data_id,
00134                         dof_idx.dimension,
00135                         dof_idx.geometry_index,
00136                         os);
00137       os << di.identity
00138          << di.interp_point;
00139       for(u_int j=0; j<fun.size(); j++)
00140         os<< fun[j](i);
00141     }
00142   }
00143 
00151   template <class FUNC>
00152     void import_fe_func(FUNC& fun,
00153                         const data_id_t& data_id) {
00154     typedef FUNC fe_func_t;
00155     typedef typename fe_func_t::fe_space_t fe_space_t;
00156     typedef typename fe_space_t::dof_info_t dof_info_t;
00157     typedef RegularMesh<fe_space_t::dim,fe_space_t::dow> mesh_t;
00158     const fe_space_t& sp = fun.femSpace();
00159     const mesh_t& mesh = dynamic_cast<const mesh_t&>(sp.mesh());
00160     AFEPack::istream<> is;
00161     u_int n_dof = sp.n_dof();
00162     for (u_int i = 0;i < n_dof;++ i) {
00163       const dof_info_t& di = sp.dofInfo(i);
00164       const DOFIndex& dof_idx = sp.dofIndex(i);
00165       double D = di.interp_point.length();
00166       get_import_stream(mesh,
00167                         data_id,
00168                         dof_idx.dimension,
00169                         dof_idx.geometry_index,
00170                         is);
00171       do {
00172         dof_info_t dof_info;
00173         is >> dof_info.identity
00174            >> dof_info.interp_point 
00175            >> fun(i);
00176         if (!(dof_info.identity == di.identity)) continue;
00177         double D1 = dof_info.interp_point.length() + D;
00178         if (distance(dof_info.interp_point, 
00179                      di.interp_point) <= 1.0e-08*D1) break;
00180       } while (1);
00181     }
00182   }
00183 
00193   template <class FUNC>
00194     void import_fe_func(std::vector<FUNC> & fun,
00195                         const data_id_t& data_id) {
00196     typedef FUNC fe_func_t;
00197     typedef typename fe_func_t::fe_space_t fe_space_t;
00198     typedef typename fe_space_t::dof_info_t dof_info_t;
00199     typedef RegularMesh<fe_space_t::dim,fe_space_t::dow> mesh_t;
00200     const fe_space_t& sp = fun[0].femSpace();
00201     const mesh_t& mesh = dynamic_cast<const mesh_t&>(sp.mesh());
00202     AFEPack::istream<> is;
00203     u_int n_dof = sp.n_dof();
00204     for (u_int i = 0;i < n_dof;++ i) {
00205       const dof_info_t& di = sp.dofInfo(i);
00206       const DOFIndex& dof_idx = sp.dofIndex(i);
00207       double D = di.interp_point.length();
00208       get_import_stream(mesh,
00209                         data_id,
00210                         dof_idx.dimension,
00211                         dof_idx.geometry_index,
00212                         is);
00213       do {
00214         dof_info_t dof_info;
00215         is >> dof_info.identity
00216            >> dof_info.interp_point;
00217         for(u_int j=0; j<fun.size(); j++)
00218           is >> fun[j](i);
00219         if (!(dof_info.identity == di.identity)) continue;
00220         double D1 = dof_info.interp_point.length() + D;
00221         if (distance(dof_info.interp_point, 
00222                      di.interp_point) <= 1.0e-08*D1) break;
00223       } while (1);
00224     }
00225   }
00226 
00227 
00228 }
00229 
00230 #endif // __MPI_Migration_h__
00231