Chaste  Release::3.4
AbstractContinuumMechanicsSolver.hpp
1 /*
2 
3 Copyright (c) 2005-2016, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34 */
35 
36 #ifndef ABSTRACTCONTINUUMMECHANICSSOLVER_HPP_
37 #define ABSTRACTCONTINUUMMECHANICSSOLVER_HPP_
38 
39 #include "ContinuumMechanicsProblemDefinition.hpp"
40 #include "CompressibilityType.hpp"
41 #include "LinearSystem.hpp"
42 #include "OutputFileHandler.hpp"
43 #include "ReplicatableVector.hpp"
44 #include "AbstractTetrahedralMesh.hpp"
45 #include "QuadraticMesh.hpp"
46 #include "DistributedQuadraticMesh.hpp"
47 #include "Warnings.hpp"
48 #include "PetscException.hpp"
49 #include "GaussianQuadratureRule.hpp"
50 #include "PetscTools.hpp"
51 #include "MechanicsEventHandler.hpp"
52 #include "CommandLineArguments.hpp"
53 #include "VtkMeshWriter.hpp"
54 
55 
61 typedef enum _ApplyDirichletBcsType
62 {
63  LINEAR_PROBLEM,
64  NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY,
65  NONLINEAR_PROBLEM_APPLY_TO_EVERYTHING
66 } ApplyDirichletBcsType;
67 
68 //forward declarations
69 template<unsigned DIM> class StressRecoveror;
70 template<unsigned DIM> class VtkNonlinearElasticitySolutionWriter;
71 
86 template<unsigned DIM>
88 {
89  friend class StressRecoveror<DIM>;
90  friend class VtkNonlinearElasticitySolutionWriter<DIM>;
91 
92 protected:
98 
101 
104 
106  std::string mOutputDirectory;
107 
110 
114  std::vector<c_vector<double,DIM> > mSpatialSolution;
115 
117  std::vector<double> mPressureSolution;
118 
119 
126  std::vector<double> mCurrentSolution;
127 
130 
133 
138  CompressibilityType mCompressibilityType;
139 
149 
153  unsigned mNumDofs;
154 
160  bool mVerbose;
161 
183 
188 
193 
198 
203 
207  void AllocateMatrixMemory();
208 
209 
235  void ApplyDirichletBoundaryConditions(ApplyDirichletBcsType type, bool symmetricProblem);
236 
269  void AddIdentityBlockForDummyPressureVariables(ApplyDirichletBcsType type);
270 
271 
289 
290 public:
299  ContinuumMechanicsProblemDefinition<DIM>& rProblemDefinition,
300  std::string outputDirectory,
301  CompressibilityType compressibilityType);
302 
305 
317  void WriteCurrentSpatialSolution(std::string fileName, std::string fileExtension, int counterToAppend=-1);
318 
331  void WriteCurrentPressureSolution( int counterToAppend=-1);
332 
338  void SetWriteOutput(bool writeOutput=true);
339 
347  void CreateVtkOutput(std::string spatialSolutionName="Spatial solution");
348 
353  std::vector<double>& rGetCurrentSolution()
354  {
355  return mCurrentSolution;
356  }
357 
362  virtual std::vector<c_vector<double,DIM> >& rGetSpatialSolution()=0;
363 
372  std::vector<double>& rGetPressures();
373 
374 };
375 
376 
377 template<unsigned DIM>
379  ContinuumMechanicsProblemDefinition<DIM>& rProblemDefinition,
380  std::string outputDirectory,
381  CompressibilityType compressibilityType)
382  : mrQuadMesh(rQuadMesh),
383  mrProblemDefinition(rProblemDefinition),
384  mOutputDirectory(outputDirectory),
385  mpOutputFileHandler(NULL),
386  mpQuadratureRule(NULL),
387  mpBoundaryQuadratureRule(NULL),
388  mCompressibilityType(compressibilityType),
389  mResidualVector(NULL),
390  mSystemLhsMatrix(NULL),
391  mPreconditionMatrix(NULL)
392 {
393  assert(DIM==2 || DIM==3);
394 
395  //Check that the mesh is Quadratic
396  QuadraticMesh<DIM>* p_quad_mesh = dynamic_cast<QuadraticMesh<DIM>* >(&rQuadMesh);
397  DistributedQuadraticMesh<DIM>* p_distributed_quad_mesh = dynamic_cast<DistributedQuadraticMesh<DIM>* >(&rQuadMesh);
398 
399  if(p_quad_mesh == NULL && p_distributed_quad_mesh == NULL)
400  {
401  EXCEPTION("Continuum mechanics solvers require a quadratic mesh");
402  }
403 
404 
405  mVerbose = (mrProblemDefinition.GetVerboseDuringSolve() ||
406  CommandLineArguments::Instance()->OptionExists("-mech_verbose") ||
407  CommandLineArguments::Instance()->OptionExists("-mech_very_verbose") );
408 
409  mWriteOutput = (mOutputDirectory != "");
410  if (mWriteOutput)
411  {
413  }
414 
415  // see dox for mProblemDimension
416  mProblemDimension = mCompressibilityType==COMPRESSIBLE ? DIM : DIM+1;
418 
420 
421  // In general the Jacobian for a mechanics problem is non-polynomial.
422  // We therefore use the highest order integration rule available.
424  // The boundary forcing terms (or tractions) are also non-polynomial in general.
425  // Again, we use the highest order integration rule available.
427 
428  mCurrentSolution.resize(mNumDofs, 0.0);
429 }
430 
431 
432 template<unsigned DIM>
434 {
435  if (mpOutputFileHandler)
436  {
437  delete mpOutputFileHandler;
438  }
439 
440  if (mpQuadratureRule)
441  {
442  delete mpQuadratureRule;
443  delete mpBoundaryQuadratureRule;
444  }
445 
446  if (mResidualVector)
447  {
448  PetscTools::Destroy(mResidualVector);
449  PetscTools::Destroy(mLinearSystemRhsVector);
450  PetscTools::Destroy(mSystemLhsMatrix);
451  PetscTools::Destroy(mPreconditionMatrix);
452  }
453 
454  if (mDirichletBoundaryConditionsVector)
455  {
456  PetscTools::Destroy(mDirichletBoundaryConditionsVector);
457  }
458 }
459 
460 template<unsigned DIM>
462  std::string fileExtension,
463  int counterToAppend)
464 {
465  // Only write output if the flag mWriteOutput has been set
466  if (!mWriteOutput)
467  {
468  return;
469  }
470 
471  if (PetscTools::AmMaster())
472  {
473  std::stringstream file_name;
474 
475  file_name << fileName;
476  if (counterToAppend >= 0)
477  {
478  file_name << "_" << counterToAppend;
479  }
480  file_name << "." << fileExtension;
481 
482  out_stream p_file = mpOutputFileHandler->OpenOutputFile(file_name.str());
483 
484  std::vector<c_vector<double,DIM> >& r_spatial_solution = rGetSpatialSolution();
485  for (unsigned i=0; i<r_spatial_solution.size(); i++)
486  {
487  // for (unsigned j=0; j<DIM; j++)
488  // {
489  // *p_file << mrQuadMesh.GetNode(i)->rGetLocation()[j] << " ";
490  // }
491  for (unsigned j=0; j<DIM; j++)
492  {
493  *p_file << r_spatial_solution[i](j) << " ";
494  }
495  *p_file << "\n";
496  }
497  p_file->close();
498  }
499  PetscTools::Barrier("WriteSpatial");
500 }
501 
502 template<unsigned DIM>
504 {
505  // Only write output if the flag mWriteOutput has been set
506  if (!mWriteOutput)
507  {
508  return;
509  }
510 
511  if (PetscTools::AmMaster())
512  {
513  std::stringstream file_name;
514 
515  file_name << "pressure";
516  if (counterToAppend >= 0)
517  {
518  file_name << "_" << counterToAppend;
519  }
520  file_name << ".txt";
521 
522  out_stream p_file = mpOutputFileHandler->OpenOutputFile(file_name.str());
523 
524  std::vector<double>& r_pressure = rGetPressures();
525  for (unsigned i=0; i<r_pressure.size(); i++)
526  {
527  for (unsigned j=0; j<DIM; j++)
528  {
529  *p_file << mrQuadMesh.GetNode(i)->rGetLocation()[j] << " ";
530  }
531 
532  *p_file << r_pressure[i] << "\n";
533  }
534  p_file->close();
535  }
536  PetscTools::Barrier("WritePressure");
537 }
538 
539 
540 template<unsigned DIM>
542 {
543  if (writeOutput && (mOutputDirectory==""))
544  {
545  EXCEPTION("Can't write output if no output directory was given in constructor");
546  }
547  mWriteOutput = writeOutput;
548 }
549 
550 // ** TO BE DEPRECATED - see #2321 **
551 template<unsigned DIM>
552 void AbstractContinuumMechanicsSolver<DIM>::CreateVtkOutput(std::string spatialSolutionName)
553 {
554  if (this->mOutputDirectory=="")
555  {
556  EXCEPTION("No output directory was given so no output was written, cannot convert to VTK format");
557  }
558 #ifdef CHASTE_VTK
559  VtkMeshWriter<DIM, DIM> mesh_writer(this->mOutputDirectory + "/vtk", "solution", true);
560 
561  mesh_writer.AddPointData(spatialSolutionName, this->rGetSpatialSolution());
562 
563  if (mCompressibilityType==INCOMPRESSIBLE)
564  {
565  mesh_writer.AddPointData("Pressure", rGetPressures());
566  }
567 
568  //Output the element attribute as cell data.
569  std::vector<double> element_attribute;
570  for(typename QuadraticMesh<DIM>::ElementIterator iter = this->mrQuadMesh.GetElementIteratorBegin();
571  iter != this->mrQuadMesh.GetElementIteratorEnd();
572  ++iter)
573  {
574  element_attribute.push_back(iter->GetAttribute());
575  }
576  mesh_writer.AddCellData("Attribute", element_attribute);
577 
578  mesh_writer.WriteFilesUsingMesh(this->mrQuadMesh);
579 #endif
580 }
581 
582 
583 template<unsigned DIM>
585 {
586  assert(mProblemDimension==DIM+1);
587 
588  mPressureSolution.clear();
589  mPressureSolution.resize(mrQuadMesh.GetNumNodes());
590 
591  for (unsigned i=0; i<mrQuadMesh.GetNumNodes(); i++)
592  {
593  mPressureSolution[i] = mCurrentSolution[mProblemDimension*i + DIM];
594  }
595  return mPressureSolution;
596 }
597 
598 
599 template<unsigned DIM>
601 {
602  assert(mProblemDimension==DIM+1);
603 
604  // For quadratic triangles, node 3 is between nodes 1 and 2, node 4 is between 0 and 2, etc
605  unsigned internal_nodes_2d[3] = {3,4,5};
606  unsigned neighbouring_vertices_2d[3][2] = { {1,2}, {2,0}, {0,1} };
607 
608  // ordering for quadratic tetrahedra
609  unsigned internal_nodes_3d[6] = {4,5,6,7,8,9};
610  unsigned neighbouring_vertices_3d[6][2] = { {0,1}, {1,2}, {0,2}, {0,3}, {1,3}, {2,3} };
611 
612  unsigned num_internal_nodes_per_element = DIM==2 ? 3 : 6;
613 
614  // loop over elements, then loop over edges.
616  iter != mrQuadMesh.GetElementIteratorEnd();
617  ++iter)
618  {
619  for(unsigned i=0; i<num_internal_nodes_per_element; i++)
620  {
621  unsigned global_index;
622  double left_val;
623  double right_val;
624 
625  if(DIM==2)
626  {
627  global_index = iter->GetNodeGlobalIndex( internal_nodes_2d[i] );
628  unsigned vertex_0_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_2d[i][0] );
629  unsigned vertex_1_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_2d[i][1] );
630  left_val = mCurrentSolution[mProblemDimension*vertex_0_global_index + DIM];
631  right_val = mCurrentSolution[mProblemDimension*vertex_1_global_index + DIM];
632  }
633  else
634  {
635  global_index = iter->GetNodeGlobalIndex( internal_nodes_3d[i] );
636  unsigned vertex_0_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_3d[i][0] );
637  unsigned vertex_1_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_3d[i][1] );
638  left_val = mCurrentSolution[mProblemDimension*vertex_0_global_index + DIM];
639  right_val = mCurrentSolution[mProblemDimension*vertex_1_global_index + DIM];
640  }
641 
642  // this line assumes the internal node is midway between the two vertices
643  mCurrentSolution[mProblemDimension*global_index + DIM] = 0.5 * (left_val + right_val);
644  }
645  }
646 }
647 
648 /*
649  * This method applies the appropriate BCs to the residual and/or linear system
650  *
651  * For the latter, and second input parameter==true, the BCs are imposed in such a way as
652  * to ensure that a symmetric linear system remains symmetric. For each row with boundary condition
653  * applied, both the row and column are zero'd and the RHS vector modified to take into account the
654  * zero'd column.
655  *
656  * Suppose we have a matrix
657  * [a b c] [x] = [ b1 ]
658  * [d e f] [y] [ b2 ]
659  * [g h i] [z] [ b3 ]
660  * and we want to apply the boundary condition x=v without losing symmetry if the matrix is
661  * symmetric. We apply the boundary condition
662  * [1 0 0] [x] = [ v ]
663  * [d e f] [y] [ b2 ]
664  * [g h i] [z] [ b3 ]
665  * and then zero the column as well, adding a term to the RHS to take account for the
666  * zero-matrix components
667  * [1 0 0] [x] = [ v ] - v[ 0 ]
668  * [0 e f] [y] [ b2 ] [ d ]
669  * [0 h i] [z] [ b3 ] [ g ]
670  * Note the last term is the first column of the matrix, with one component zeroed, and
671  * multiplied by the boundary condition value. This last term is then stored in
672  * mDirichletBoundaryConditionsVector, and in general is equal to:
673  * SUM_{d=1..D} v_d a'_d
674  * where v_d is the boundary value of boundary condition d (d an index into the matrix),
675  * and a'_d is the dth-column of the matrix but with the d-th component zeroed, and where
676  * there are D boundary conditions
677  */
678 template<unsigned DIM>
679 void AbstractContinuumMechanicsSolver<DIM>::ApplyDirichletBoundaryConditions(ApplyDirichletBcsType type, bool symmetricProblem)
680 {
681  std::vector<unsigned> rows;
682  std::vector<double> values;
683 
684  // Whether to apply symmetrically, ie alter columns as well as rows (see comment above)
685  bool applySymmetrically = (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY) && symmetricProblem;
686 
687  if (applySymmetrically)
688  {
689  if(mDirichletBoundaryConditionsVector==NULL)
690  {
691  VecDuplicate(mResidualVector, &mDirichletBoundaryConditionsVector);
692  }
693 
694  PetscVecTools::Zero(mDirichletBoundaryConditionsVector);
695  PetscMatTools::Finalise(mSystemLhsMatrix);
696  }
697 
699  // collect the entries to be altered
701 
702  for (unsigned i=0; i<mrProblemDefinition.rGetDirichletNodes().size(); i++)
703  {
704  unsigned node_index = mrProblemDefinition.rGetDirichletNodes()[i];
705 
706  for (unsigned j=0; j<DIM; j++)
707  {
708  double dirichlet_val = mrProblemDefinition.rGetDirichletNodeValues()[i](j);
709 
711  {
712  double val;
713  unsigned dof_index = mProblemDimension*node_index+j;
714 
715  if(type == LINEAR_PROBLEM)
716  {
717  val = dirichlet_val;
718  }
719  else
720  {
721  // The boundary conditions on the NONLINEAR SYSTEM are x=boundary_values
722  // on the boundary nodes. However:
723  // The boundary conditions on the LINEAR SYSTEM at each Newton step (Ju=f,
724  // where J is the Jacobian, u the negative update vector and f is the residual) is:
725  // u=current_soln-boundary_values on the boundary nodes
726  val = mCurrentSolution[dof_index] - dirichlet_val;
727  }
728  rows.push_back(dof_index);
729  values.push_back(val);
730  }
731  }
732  }
733 
735  // do the alterations
737 
738  if (applySymmetrically)
739  {
740  // Modify the matrix columns
741  for (unsigned i=0; i<rows.size(); i++)
742  {
743  unsigned col = rows[i];
744  double minus_value = -values[i];
745 
746  // Get a vector which will store the column of the matrix (column d, where d is
747  // the index of the row (and column) to be altered for the boundary condition.
748  // Since the matrix is symmetric when get row number "col" and treat it as a column.
749  // PETSc uses compressed row format and therefore getting rows is far more efficient
750  // than getting columns.
751  Vec matrix_col = PetscMatTools::GetMatrixRowDistributed(mSystemLhsMatrix,col);
752 
753  // Zero the correct entry of the column
754  PetscVecTools::SetElement(matrix_col, col, 0.0);
755 
756  // Set up the RHS Dirichlet boundary conditions vector
757  // E.g. for a boundary node at the zeroth node (x_0 = value), this is equal to
758  // -value*[0 a_21 a_31 .. a_N1]
759  // and will be added to the RHS.
760  PetscVecTools::AddScaledVector(mDirichletBoundaryConditionsVector, matrix_col, minus_value);
761  PetscTools::Destroy(matrix_col);
762  }
763  }
764 
765  if (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY) // i.e. doing a whole linear system
766  {
767  // Now zero the appropriate rows and columns of the matrix. If the matrix is symmetric we apply the
768  // boundary conditions in a way the symmetry isn't lost (rows and columns). If not only the row is
769  // zeroed.
770  if (applySymmetrically)
771  {
772  PetscMatTools::ZeroRowsAndColumnsWithValueOnDiagonal(mSystemLhsMatrix, rows, 1.0);
773  PetscMatTools::ZeroRowsAndColumnsWithValueOnDiagonal(mPreconditionMatrix, rows, 1.0);
774 
775  // Apply the RHS boundary conditions modification if required.
776  PetscVecTools::AddScaledVector(mLinearSystemRhsVector, mDirichletBoundaryConditionsVector, 1.0);
777  }
778  else
779  {
780  PetscMatTools::ZeroRowsWithValueOnDiagonal(mSystemLhsMatrix, rows, 1.0);
781  PetscMatTools::ZeroRowsWithValueOnDiagonal(mPreconditionMatrix, rows, 1.0);
782  }
783  }
784 
785  if (type!=LINEAR_PROBLEM)
786  {
787  for (unsigned i=0; i<rows.size(); i++)
788  {
789  PetscVecTools::SetElement(mResidualVector, rows[i], values[i]);
790  }
791  }
792 
793  if (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY)
794  {
795  for (unsigned i=0; i<rows.size(); i++)
796  {
797  PetscVecTools::SetElement(mLinearSystemRhsVector, rows[i], values[i]);
798  }
799  }
800 }
801 
802 
803 template<unsigned DIM>
805 {
806  assert(mCompressibilityType==INCOMPRESSIBLE);
807 
808  int lo, hi;
809  VecGetOwnershipRange(mResidualVector, &lo, &hi);
810 
811  for(unsigned i=0; i<mrQuadMesh.GetNumNodes(); i++)
812  {
813  if(mrQuadMesh.GetNode(i)->IsInternal())
814  {
815  unsigned row = (DIM+1)*i + DIM; // DIM+1 is the problem dimension
816  if(lo <= (int)row && (int)row < hi)
817  {
818  if(type!=LINEAR_PROBLEM)
819  {
820  PetscVecTools::SetElement(mResidualVector, row, mCurrentSolution[row]-0.0);
821  }
822  if(type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY) // ie doing a whole linear system
823  {
824  double rhs_vector_val = type==LINEAR_PROBLEM ? 0.0 : mCurrentSolution[row]-0.0;
825  PetscVecTools::SetElement(mLinearSystemRhsVector, row, rhs_vector_val);
826  // this assumes the row is already zero, which is should be..
827  PetscMatTools::SetElement(mSystemLhsMatrix, row, row, 1.0);
828  PetscMatTools::SetElement(mPreconditionMatrix, row, row, 1.0);
829  }
830  }
831  }
832  }
833 }
834 
835 
836 
837 template<unsigned DIM>
839 {
840  Vec template_vec = mrQuadMesh.GetDistributedVectorFactory()->CreateVec(mProblemDimension);
841 
843  // three vectors
845  VecDuplicate(template_vec, &mResidualVector);
846  VecDuplicate(mResidualVector, &mLinearSystemRhsVector);
847  // the one is only allocated if it will be needed (in ApplyDirichletBoundaryConditions),
848  // depending on whether the matrix is kept symmetric.
849  mDirichletBoundaryConditionsVector = NULL;
850  PetscTools::Destroy(template_vec);
851 
853  // two matrices
855 
856  int lo, hi;
857  VecGetOwnershipRange(mResidualVector, &lo, &hi);
858  PetscInt local_size = hi - lo;
859 
860 
861  if (DIM==2)
862  {
863  // 2D: N elements around a point => 7N+3 non-zeros in that row? Assume N<=10 (structured mesh would have N_max=6) => 73.
864  unsigned num_non_zeros = std::min(75u, mNumDofs);
865 
866  PetscTools::SetupMat(mSystemLhsMatrix, mNumDofs, mNumDofs, num_non_zeros, local_size, local_size);
867  PetscTools::SetupMat(mPreconditionMatrix, mNumDofs, mNumDofs, num_non_zeros, local_size, local_size);
868  }
869  else
870  {
871  assert(DIM==3);
872 
873  // in 3d we get the number of containing elements for each node and use that to obtain an upper bound
874  // for the number of non-zeros for each DOF associated with that node.
875 
876  int* num_non_zeros_each_row = new int[mNumDofs];
877  for (unsigned i=0; i<mNumDofs; i++)
878  {
879  num_non_zeros_each_row[i] = 0;
880  }
881 
882  for (typename AbstractMesh<DIM,DIM>::NodeIterator iter = mrQuadMesh.GetNodeIteratorBegin();
883  iter != mrQuadMesh.GetNodeIteratorEnd();
884  ++iter)
885  {
886  // this upper bound neglects the fact that two containing elements will share the same nodes..
887  // 4 = max num dofs associated with this node
888  // 30 = 3*9+3 = 3 dimensions x 9 other nodes on this element + 3 vertices with a pressure unknown
889  unsigned num_non_zeros_upper_bound = 4 + 30*iter->GetNumContainingElements();
890 
891  num_non_zeros_upper_bound = std::min(num_non_zeros_upper_bound, mNumDofs);
892 
893  unsigned i = iter->GetIndex();
894 
895  num_non_zeros_each_row[mProblemDimension*i + 0] = num_non_zeros_upper_bound;
896  num_non_zeros_each_row[mProblemDimension*i + 1] = num_non_zeros_upper_bound;
897  num_non_zeros_each_row[mProblemDimension*i + 2] = num_non_zeros_upper_bound;
898 
899  if (mCompressibilityType==INCOMPRESSIBLE)
900  {
901  if(!iter->IsInternal())
902  {
903  num_non_zeros_each_row[mProblemDimension*i + 3] = num_non_zeros_upper_bound;
904  }
905  else
906  {
907  num_non_zeros_each_row[mProblemDimension*i + 3] = 1;
908  }
909  }
910  }
911 
912  // NOTE: PetscTools::SetupMat() or the below creates a MATAIJ matrix, which means the matrix will
913  // be of type MATSEQAIJ if num_procs=1 and MATMPIAIJ otherwise. In the former case
914  // MatSeqAIJSetPreallocation MUST be called [MatMPIAIJSetPreallocation will have
915  // no effect (silently)], and vice versa in the latter case
916 
919  //PetscTools::SetupMat(mSystemLhsMatrix, mNumDofs, mNumDofs, 0, PETSC_DECIDE, PETSC_DECIDE);
920  //PetscTools::SetupMat(mPreconditionMatrix, mNumDofs, mNumDofs, 0, PETSC_DECIDE, PETSC_DECIDE);
922 
923  // possible todo: create a PetscTools::SetupMatNoAllocation()
924 
925 
926 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
927  MatCreate(PETSC_COMM_WORLD,local_size,local_size,mNumDofs,mNumDofs,&mSystemLhsMatrix);
928  MatCreate(PETSC_COMM_WORLD,local_size,local_size,mNumDofs,mNumDofs,&mPreconditionMatrix);
929 #else //New API
930  MatCreate(PETSC_COMM_WORLD,&mSystemLhsMatrix);
931  MatCreate(PETSC_COMM_WORLD,&mPreconditionMatrix);
932  MatSetSizes(mSystemLhsMatrix,local_size,local_size,mNumDofs,mNumDofs);
933  MatSetSizes(mPreconditionMatrix,local_size,local_size,mNumDofs,mNumDofs);
934 #endif
935 
937  {
938  MatSetType(mSystemLhsMatrix, MATSEQAIJ);
939  MatSetType(mPreconditionMatrix, MATSEQAIJ);
940  MatSeqAIJSetPreallocation(mSystemLhsMatrix, PETSC_DEFAULT, num_non_zeros_each_row);
941  MatSeqAIJSetPreallocation(mPreconditionMatrix, PETSC_DEFAULT, num_non_zeros_each_row);
942  }
943  else
944  {
945  int* num_non_zeros_each_row_in_diag = new int[local_size];
946  int* num_non_zeros_each_row_off_diag = new int[local_size];
947  for (unsigned i=0; i<unsigned(local_size); i++)
948  {
949  num_non_zeros_each_row_in_diag[i] = num_non_zeros_each_row[lo+i];
950  num_non_zeros_each_row_off_diag[i] = num_non_zeros_each_row[lo+i];
951  // In the on process ("diagonal block") there cannot be more non-zero columns specified than there are rows
952  if(num_non_zeros_each_row_in_diag[i] > local_size)
953  {
954  num_non_zeros_each_row_in_diag[i] = local_size;
955  }
956  }
957 
958  MatSetType(mSystemLhsMatrix, MATMPIAIJ);
959  MatSetType(mPreconditionMatrix, MATMPIAIJ);
960  MatMPIAIJSetPreallocation(mSystemLhsMatrix, PETSC_DEFAULT, num_non_zeros_each_row_in_diag, PETSC_DEFAULT, num_non_zeros_each_row_off_diag);
961  MatMPIAIJSetPreallocation(mPreconditionMatrix, PETSC_DEFAULT, num_non_zeros_each_row_in_diag, PETSC_DEFAULT, num_non_zeros_each_row_off_diag);
962  }
963 
964  MatSetFromOptions(mSystemLhsMatrix);
965  MatSetFromOptions(mPreconditionMatrix);
966 #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
967  MatSetOption(mSystemLhsMatrix, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE);
968  MatSetOption(mPreconditionMatrix, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE);
969 #else
970  MatSetOption(mSystemLhsMatrix, MAT_IGNORE_OFF_PROC_ENTRIES);
971  MatSetOption(mPreconditionMatrix, MAT_IGNORE_OFF_PROC_ENTRIES);
972 #endif
973 
974  //unsigned total_non_zeros = 0;
975  //for (unsigned i=0; i<mNumDofs; i++)
976  //{
977  // total_non_zeros += num_non_zeros_each_row[i];
978  //}
979  //std::cout << total_non_zeros << " versus " << 500*mNumDofs << "\n" << std::flush;
980 
981  delete [] num_non_zeros_each_row;
982  }
983 }
984 #endif // ABSTRACTCONTINUUMMECHANICSSOLVER_HPP_
ElementIterator GetElementIteratorBegin(bool skipDeletedElements=true)
void ApplyDirichletBoundaryConditions(ApplyDirichletBcsType type, bool symmetricProblem)
GaussianQuadratureRule< DIM > * mpQuadratureRule
void WriteCurrentPressureSolution(int counterToAppend=-1)
AbstractTetrahedralMesh< DIM, DIM > & mrQuadMesh
static void Barrier(const std::string callerId="")
Definition: PetscTools.cpp:134
virtual std::vector< c_vector< double, DIM > > & rGetSpatialSolution()=0
void WriteCurrentSpatialSolution(std::string fileName, std::string fileExtension, int counterToAppend=-1)
#define EXCEPTION(message)
Definition: Exception.hpp:143
static void SetElement(Vec vector, PetscInt row, double value)
static bool AmMaster()
Definition: PetscTools.cpp:120
void AddCellData(std::string name, std::vector< double > data)
std::vector< c_vector< double, DIM > > mSpatialSolution
virtual unsigned GetNumNodes() const
AbstractContinuumMechanicsSolver(AbstractTetrahedralMesh< DIM, DIM > &rQuadMesh, ContinuumMechanicsProblemDefinition< DIM > &rProblemDefinition, std::string outputDirectory, CompressibilityType compressibilityType)
void AddIdentityBlockForDummyPressureVariables(ApplyDirichletBcsType type)
bool OptionExists(const std::string &rOption)
void WriteFilesUsingMesh(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, bool keepOriginalElementIndexing=true)
static Vec GetMatrixRowDistributed(Mat matrix, unsigned rowIndex)
static bool IsSequential()
Definition: PetscTools.cpp:91
ContinuumMechanicsProblemDefinition< DIM > & mrProblemDefinition
static void SetElement(Mat matrix, PetscInt row, PetscInt col, double value)
static void ZeroRowsWithValueOnDiagonal(Mat matrix, std::vector< unsigned > &rRows, double diagonalValue)
static void ZeroRowsAndColumnsWithValueOnDiagonal(Mat matrix, std::vector< unsigned > rowColIndices, double diagonalValue)
static void AddScaledVector(Vec y, Vec x, double scaleFactor)
static void SetupMat(Mat &rMat, int numRows, int numColumns, unsigned rowPreallocation, int numLocalRows=PETSC_DECIDE, int numLocalColumns=PETSC_DECIDE, bool ignoreOffProcEntries=true, bool newAllocationError=true)
Definition: PetscTools.cpp:268
NodeIterator GetNodeIteratorBegin(bool skipDeletedNodes=true)
static void Zero(Vec vector)
void CreateVtkOutput(std::string spatialSolutionName="Spatial solution")
static void Destroy(Vec &rVec)
Definition: PetscTools.hpp:351
static void Finalise(Mat matrix)
static CommandLineArguments * Instance()
GaussianQuadratureRule< DIM-1 > * mpBoundaryQuadratureRule
void AddPointData(std::string name, std::vector< double > data)