• Main Page
  • Related Pages
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

Interface.hpp

Go to the documentation of this file.
00001 /* -*- mode: c++ -*-
00002 
00003   This file is part of the LifeV Applications.
00004 
00005   Author(s):
00006        Date: 2009-03-24
00007 
00008   Copyright (C) 2009 EPFL
00009 
00010   This program is free software; you can redistribute it and/or modify
00011   it under the terms of the GNU General Public License as published by
00012   the Free Software Foundation; either version 2.1 of the License, or
00013   (at your option) any later version.
00014 
00015   This program is distributed in the hope that it will be useful, but
00016   WITHOUT ANY WARRANTY; without even the implied warranty of
00017   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018   General Public License for more details.
00019 
00020   You should have received a copy of the GNU General Public License
00021   along with this program; if not, write to the Free Software
00022   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023   USA
00024 */
00025 
00026 
00027 // ===================================================
00029 // ===================================================
00030 
00031 //#include <boost/program_options.hpp>
00032 
00033 #define velocity_solver_init_mpi velocity_solver_init_mpi_
00034 #define velocity_solver_finalize velocity_solver_finalize_
00035 #define velocity_solver_init_l1l2 velocity_solver_init_l1l2_
00036 #define velocity_solver_solve_l1l2 velocity_solver_solve_l1l2_
00037 #define velocity_solver_init_fo velocity_solver_init_fo_
00038 #define velocity_solver_solve_fo velocity_solver_solve_fo_
00039 #define velocity_solver_init_stokes velocity_solver_init_stokes_
00040 #define velocity_solver_solve_stokes velocity_solver_solve_stokes_
00041 #define velocity_solver_compute_2d_grid velocity_solver_compute_2d_grid_
00042 #define velocity_solver_set_grid_data velocity_solver_set_grid_data_
00043 #define velocity_solver_extrude_3d_grid velocity_solver_extrude_3d_grid_
00044 #define velocity_solver_export_l1l2_velocity velocity_solver_export_l1l2_velocity_ 
00045 #define velocity_solver_export_2d_data velocity_solver_export_2d_data_ 
00046 #define velocity_solver_export_fo_velocity velocity_solver_export_fo_velocity_
00047 #define velocity_solver_estimate_SS_SMB velocity_solver_estimate_ss_smb_
00048 /*
00049 #include "Extrude3DMesh.hpp"
00050 /*/
00051 #include <vector>
00052 #include <mpi.h>
00053 #include <list>
00054 #include <iostream>
00055 #include <limits>
00056 #include <cmath>
00057 
00058 //enum ordering{LayerWise, ColumnWise};
00059 
00060 typedef unsigned int ID;
00061 typedef unsigned int UInt;
00062 const ID NotAnId = std::numeric_limits<int>::max();
00063 //*/
00064 // ===================================================
00066 // ===================================================
00067 extern "C" {
00068 
00069 // 1
00070 int velocity_solver_init_mpi(int *fComm);
00071 
00072 void velocity_solver_finalize();
00073 
00074 void velocity_solver_init_l1l2(double const * levelsRatio);
00075 
00076 // 5
00077 void velocity_solver_init_fo(double const * levelsRatio);
00078 
00079 void velocity_solver_solve_l1l2(double const * lowerSurface_F, double const * thickness_F,
00080                double const * beta_F, double const * temperature_F,
00081                double * u_normal_F = 0,
00082                double * heatIntegral_F = 0 , double * viscosity_F = 0);
00083 
00084 // 6
00085 void velocity_solver_solve_fo(double const * lowerSurface_F, double const * thickness_F,
00086                            double const * beta_F, double const * temperature_F,
00087                            double * u_normal_F = 0,
00088                            double * heatIntegral_F = 0 , double * viscosity_F = 0);
00089 
00090 
00091 // 3
00092 void velocity_solver_compute_2d_grid(int const * verticesMask_F);
00093 
00094 
00095 // 2
00096 void velocity_solver_set_grid_data(int const * _nCells_F, int const * _nEdges_F, int const * _nVertices_F, int const * _nLayers,
00097                                  int const * _nCellsSolve_F, int const * _nEdgesSolve_F, int const * _nVerticesSolve_F, int const* _maxNEdgesOnCell_F,
00098                                  double const * radius_F,
00099                                  int const * _cellsOnEdge_F, int const * _cellsOnVertex_F, int const * _verticesOnCell_F, int const * _verticesOnEdge_F, int const * _edgesOnCell_F,
00100                                  int const* _nEdgesOnCells_F, int const * _indexToCellID_F,
00101                                  double const *  _xCell_F, double const *  _yCell_F, double const *  _zCell_F, double const *  _areaTriangle_F,
00102                                  int const * sendCellsArray_F, int const * recvCellsArray_F,
00103                                  int const * sendEdgesArray_F, int const * recvEdgesArray_F,
00104                                  int const * sendVerticesArray_F, int const * recvVerticesArray_F);
00105 
00106 // 4
00107 void velocity_solver_extrude_3d_grid(double const * levelsRatio_F, double const * lowerSurface_F, double const * thickness_F);
00108 
00109 void velocity_solver_export_l1l2_velocity();
00110 
00111 void velocity_solver_export_fo_velocity();
00112 
00113 
00114 }
00115 
00116 struct exchange{
00117             const int procID;
00118             const std::vector<int> vec;
00119             mutable std::vector<int> buffer;
00120       mutable std::vector<double> doubleBuffer;
00121             mutable MPI_Request reqID;
00122 
00123             exchange(int _procID, int const *  vec_first, int const *  vec_last, int fieldDim=1);
00124         };
00125 
00126 typedef std::list<exchange> exchangeList_Type;
00127 
00128 exchangeList_Type unpackMpiArray(int const * array);
00129 
00130 bool isGhostTriangle(int i, double relTol = 1e-1);
00131 
00132 double signedTriangleArea(const double* x, const double* y);
00133 
00134 double signedTriangleArea(const double* x, const double* y, const double* z);
00135 
00136 void import2DFields(double const * lowerSurface_F, double const * thickness_F,
00137                       double const * beta_F=0, double eps=0);
00138 
00139 std::vector<int> extendMaskByOneLayer(int const* verticesMask_F);
00140 
00141 void extendMaskByOneLayer(int const* verticesMask_F, std::vector<int>& extendedFVerticesMask);
00142 
00143 void importP0Temperature(double const * temperature_F);
00144 
00145 void get_tetraP1_velocity_on_FEdges(double * uNormal, const std::vector<double>& velocityOnVertices, const std::vector<int>& edgeToFEdge, const std::vector<int>& mpasIndexToVertexID);
00146 
00147 void get_prism_velocity_on_FEdges(double * uNormal, const std::vector<double>& velocityOnVertices, const std::vector<int>& edgeToFEdge);
00148 
00149 void createReverseCellsExchangeLists(exchangeList_Type& sendListReverse_F, exchangeList_Type& receiveListReverse_F, const std::vector<int>& fVertexToTriangleID, const std::vector<int>& fCellToVertexID);
00150 
00151 void createReverseEdgesExchangeLists(exchangeList_Type& sendListReverse_F, exchangeList_Type& receiveListReverse_F, const std::vector<int>& fVertexToTriangleID, const std::vector<int>& fEdgeToEdgeID);
00152 
00153 void mapCellsToVertices(const std::vector<double>& velocityOnCells, std::vector<double>& velocityOnVertices, int fieldDim, int numLayers, int ordering);
00154 
00155 void mapVerticesToCells(const std::vector<double>& velocityOnVertices, double* velocityOnCells, int fieldDim, int numLayers, int ordering);
00156 
00157 void createReducedMPI(int nLocalEntities, MPI_Comm& reduced_comm_id);
00158 
00159 void computeLocalOffset(int nLocalEntities, int& localOffset, int& nGlobalEntities);
00160 
00161 void getProcIds(std::vector<int>& field,int const * recvArray);
00162 
00163 void getProcIds(std::vector<int>& field, exchangeList_Type const * recvList);
00164 
00165 void allToAll(std::vector<int>& field, int const * sendArray, int const * recvArray, int fieldDim=1);
00166 
00167 void allToAll(std::vector<int>& field, exchangeList_Type const * sendList, exchangeList_Type const * recvList, int fieldDim=1);
00168 
00169 void allToAll(double* field, exchangeList_Type const * sendList, exchangeList_Type const * recvList, int fieldDim=1);
00170 
00171 
00172 int initialize_iceProblem(int nTriangles);
00173 
00174 
00175 
00176 
00177 
00178 
00179 

Generated on Wed Mar 26 2014 18:36:39 for Albany: a Trilinos-based PDE code by  doxygen 1.7.1