ImplicitCardiacMechanicsAssembler.hpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2009
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
00026 
00027 */
00028 
00029 
00030 #ifndef IMPLICITCARDIACMECHANICSASSEMBLER_HPP_
00031 #define IMPLICITCARDIACMECHANICSASSEMBLER_HPP_
00032 
00033 #include "NonlinearElasticityAssembler.hpp"
00034 #include "QuadraticBasisFunction.hpp"
00035 #include "LinearBasisFunction.hpp"
00036 #include "NhsSystemWithImplicitSolver.hpp"
00037 #include "NashHunterPoleZeroLaw.hpp"
00038 #include "LogFile.hpp"
00039 #include <cfloat>
00040 
00041 
00050 template<unsigned DIM>
00051 class ImplicitCardiacMechanicsAssembler : public NonlinearElasticityAssembler<DIM>
00052 {
00053 friend class TestImplicitCardiacMechanicsAssembler;
00054 
00055 private:
00056     static const unsigned STENCIL_SIZE = NonlinearElasticityAssembler<DIM>::STENCIL_SIZE;
00057     static const unsigned NUM_NODES_PER_ELEMENT = NonlinearElasticityAssembler<DIM>::NUM_NODES_PER_ELEMENT;
00058     static const unsigned NUM_VERTICES_PER_ELEMENT = NonlinearElasticityAssembler<DIM>::NUM_VERTICES_PER_ELEMENT;
00064     std::vector<NhsSystemWithImplicitSolver> mCellMechSystems;
00065     
00070     std::vector<double> mLambdaLastTimeStep;
00071     
00076     std::vector<double> mLambda;
00077 
00079     double mCurrentTime;
00081     double mNextTime;
00083     double mOdeTimestep;
00084 
00086     bool mAllocatedMaterialLawMemory;
00087     
00089     unsigned mTotalQuadPoints;
00090     
00091 public:
00100     ImplicitCardiacMechanicsAssembler(QuadraticMesh<DIM>* pQuadMesh,
00101                                       std::string outputDirectory,
00102                                       std::vector<unsigned>& rFixedNodes,
00103                                       AbstractIncompressibleMaterialLaw<DIM>* pMaterialLaw = NULL);
00104 
00108     ~ImplicitCardiacMechanicsAssembler();
00109     
00111     unsigned GetTotalNumQuadPoints();
00112     
00114     GaussianQuadratureRule<DIM>* GetQuadratureRule();
00115 
00122     void SetIntracellularCalciumConcentrations(std::vector<double>& caI);
00123 
00131     std::vector<double>& rGetLambda();
00132 
00140     void Solve(double currentTime, double nextTime, double odeTimestep);
00141 
00142 
00143 private:
00144 
00152     void AssembleOnElement(Element<DIM, DIM>& rElement,
00153                            c_matrix<double,STENCIL_SIZE,STENCIL_SIZE>& rAElem,
00154                            c_matrix<double,STENCIL_SIZE,STENCIL_SIZE>& rAElemPrecond,
00155                            c_vector<double,STENCIL_SIZE>& rBElem,
00156                            bool assembleResidual,
00157                            bool assembleJacobian);
00158 };
00159 
00160 
00167 //public:
00168 //    std::vector<std::vector<unsigned> > mNodesContainedInElement;
00169 //
00170 //    void ComputeElementsContainingNodes(TetrahedralMesh<DIM,DIM>* pOtherMesh)
00171 //    {
00172 //        assert(DIM==2);
00173 //
00174 //        mNodesContainedInElement.resize(this->mpMesh->n_active_cells());
00175 //
00176 //        unsigned element_number = 0;
00177 //        typename DoFHandler<DIM>::active_cell_iterator  element_iter = this->mDofHandler.begin_active();
00178 //
00179 //        while (element_iter!=this->mDofHandler.end())
00180 //        {
00181 //            double xmin = element_iter->vertex(0)(0);
00182 //            double xmax = element_iter->vertex(1)(0);
00183 //            double ymin = element_iter->vertex(0)(1);
00184 //            double ymax = element_iter->vertex(3)(1);
00185 //
00186 //            assert(element_iter->vertex(2)(0)==xmax);
00187 //            assert(element_iter->vertex(2)(1)==ymax);
00188 //
00189 //            for(unsigned i=0; i<pOtherMesh->GetNumNodes(); i++)
00190 //            {
00191 //                double x = pOtherMesh->GetNode(i)->rGetLocation()[0];
00192 //                double y = pOtherMesh->GetNode(i)->rGetLocation()[1];
00193 //                if((x>=xmin) && (x<=xmax) && (y>=ymin) && (y<=ymax))
00194 //                {
00195 //                    mNodesContainedInElement[element_number].push_back(i);
00196 //                }
00197 //            }
00198 //
00199 //            element_iter++;
00200 //            element_number++;
00201 //        }
00202 //    }
00203 //
00204 //    void WriteLambda(std::string directory, std::string fileName)
00205 //    {
00206 //        OutputFileHandler handler(directory,false);
00207 //        out_stream p_file = handler.OpenOutputFile(fileName);
00208 //
00209 //        std::vector<std::vector<double> > quad_point_posns
00210 //           = FiniteElasticityTools<DIM>::GetQuadPointPositions(*(this->mpMesh), this->GetNumQuadPointsInEachDimension());
00211 //
00212 //
00213 //        for(unsigned i=0; i<quad_point_posns.size(); i++)
00214 //        {
00215 //            (*p_file) << quad_point_posns[i][0] << " " << quad_point_posns[i][1] << " "
00216 //                      << mCellMechSystems[i].GetLambda() << "\n";
00217 //        }
00218 //    }
00219 //
00220 //
00221 //    void CalculateCinverseAtNodes(TetrahedralMesh<DIM,DIM>* pOtherMesh, std::vector<std::vector<double> >& rValuesAtNodes)
00222 //    {
00223 //        assert(DIM==2);
00224 //        rValuesAtNodes.resize(pOtherMesh->GetNumNodes());
00225 //
00226 //        unsigned element_number = 0;
00227 //
00228 //        static QTrapez<DIM>   trapezoid_quadrature_formula; //trapeziod rule - values at NODES
00229 //        const unsigned n_q_points = trapezoid_quadrature_formula.n_quadrature_points;
00230 //
00231 //        FEValues<DIM> fe_values(this->mFeSystem, trapezoid_quadrature_formula,
00232 //                                UpdateFlags(update_values    |
00233 //                                            update_gradients |
00234 //                                            update_q_points  |     // needed for interpolating u and u' on the quad point
00235 //                                            update_JxW_values));
00236 //
00237 //        std::vector< Vector<double> >                  local_solution_values(n_q_points);
00238 //        std::vector< std::vector< Tensor<1,DIM> > >    local_solution_gradients(n_q_points);
00239 //
00240 //        for (unsigned q_point=0; q_point<n_q_points; q_point++)
00241 //        {
00242 //            local_solution_values[q_point].reinit(DIM+1);
00243 //            local_solution_gradients[q_point].resize(DIM+1);
00244 //        }
00245 //
00246 //
00247 //        Tensor<2,DIM> identity;
00248 //        for (unsigned i=0; i<DIM; i++)
00249 //        {
00250 //            for (unsigned j=0; j<DIM; j++)
00251 //            {
00252 //                identity[i][j] = i==j ? 1.0 : 0.0;
00253 //            }
00254 //        }
00255 //
00256 //        typename DoFHandler<DIM>::active_cell_iterator  element_iter = this->mDofHandler.begin_active();
00257 //
00258 //        while (element_iter!=this->mDofHandler.end())
00259 //        {
00260 //            double xmin = element_iter->vertex(0)(0);
00261 //            double xmax = element_iter->vertex(1)(0);
00262 //            double ymin = element_iter->vertex(0)(1);
00263 //            double ymax = element_iter->vertex(3)(1);
00264 //            assert(element_iter->vertex(2)(0)==xmax);
00265 //            assert(element_iter->vertex(2)(1)==ymax);
00266 //
00267 //            fe_values.reinit(element_iter); // compute fe values for this element
00268 //            fe_values.get_function_values(this->mCurrentSolution, local_solution_values);
00269 //            fe_values.get_function_grads(this->mCurrentSolution, local_solution_gradients);
00270 //
00271 //            std::vector<Point<DIM> > quad_points =fe_values.get_quadrature_points();
00272 //
00273 //
00274 //            AbstractIncompressibleMaterialLaw<DIM>* p_material_law = this->GetMaterialLawForElement(element_iter);
00275 //
00276 //            std::vector<Tensor<2,DIM> > inv_C_at_nodes(4);// 4=2^DIM
00277 //
00278 //            for (unsigned q_point=0; q_point<n_q_points; q_point++)
00279 //            {
00280 //                const std::vector< Tensor<1,DIM> >& grad_u_p = local_solution_gradients[q_point];
00281 //                static Tensor<2,DIM> F;
00282 //                static Tensor<2,DIM> C;
00283 //
00284 //                for (unsigned i=0; i<DIM; i++)
00285 //                {
00286 //                    for (unsigned j=0; j<DIM; j++)
00287 //                    {
00288 //                        F[i][j] = identity[i][j] + grad_u_p[i][j];
00289 //                    }
00290 //                }
00291 //
00292 //                C = transpose(F) * F;
00293 //                inv_C_at_nodes[q_point] = invert(C);
00294 //            }
00295 //
00296 //
00306 //
00307 //
00308 //
00309 //            for(unsigned j=0; j<mNodesContainedInElement[element_number].size(); j++)
00310 //            {
00311 //                unsigned node_num = mNodesContainedInElement[element_number][j];
00312 //                double x = pOtherMesh->GetNode(node_num)->rGetLocation()[0];
00313 //                double y = pOtherMesh->GetNode(node_num)->rGetLocation()[1];
00314 //
00315 //                assert((x>=xmin) && (x<=xmax) && (y>=ymin) && (y<=ymax));
00316 //                double xi  = (x-xmin)/(xmax-xmin);
00317 //                double eta = (y-ymin)/(ymax-ymin);
00318 //                assert((0<=xi) && (x<=1) && (0<=eta) && (eta<=1));
00319 //
00320 //                rValuesAtNodes[node_num][0] = InterpolateCinverse(xi,eta,inv_C_at_nodes,0,0);
00321 //                rValuesAtNodes[node_num][1] = InterpolateCinverse(xi,eta,inv_C_at_nodes,0,1);
00322 //                rValuesAtNodes[node_num][2] = InterpolateCinverse(xi,eta,inv_C_at_nodes,1,1);
00323 //            }
00324 //
00325 //
00326 //            element_iter++;
00327 //            element_number++;
00328 //        }
00329 //    }
00330 //
00331 //
00332 //    double InterpolateCinverse(const double xi, const double eta,
00333 //                               const std::vector<Tensor<2,DIM> >& inverseCAtNodes,
00334 //                               unsigned i, unsigned j)
00335 //    {
00336 //        return    inverseCAtNodes[0][i][j] * (1-xi) * (1-eta)
00337 //                + inverseCAtNodes[1][i][j] * (1-xi) *   eta
00338 //                + inverseCAtNodes[2][i][j] *   xi   * (1-eta)
00339 //                + inverseCAtNodes[3][i][j] *   xi   *   eta;
00340 //    }
00341 
00342 
00343 #endif /*IMPLICITCARDIACMECHANICSASSEMBLER_HPP_*/

Generated on Wed Mar 18 12:51:52 2009 for Chaste by  doxygen 1.5.5