00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036 #ifndef ABSTRACTCONTINUUMMECHANICSSOLVER_HPP_
00037 #define ABSTRACTCONTINUUMMECHANICSSOLVER_HPP_
00038
00039 #include "ContinuumMechanicsProblemDefinition.hpp"
00040 #include "CompressibilityType.hpp"
00041 #include "LinearSystem.hpp"
00042 #include "OutputFileHandler.hpp"
00043 #include "ReplicatableVector.hpp"
00044 #include "AbstractTetrahedralMesh.hpp"
00045 #include "QuadraticMesh.hpp"
00046 #include "DistributedQuadraticMesh.hpp"
00047 #include "Warnings.hpp"
00048 #include "PetscException.hpp"
00049 #include "GaussianQuadratureRule.hpp"
00050 #include "PetscTools.hpp"
00051 #include "MechanicsEventHandler.hpp"
00052 #include "CommandLineArguments.hpp"
00053 #include "VtkMeshWriter.hpp"
00054
00055
00061 typedef enum _ApplyDirichletBcsType
00062 {
00063 LINEAR_PROBLEM,
00064 NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY,
00065 NONLINEAR_PROBLEM_APPLY_TO_EVERYTHING
00066 } ApplyDirichletBcsType;
00067
00068
00069 template<unsigned DIM> class StressRecoveror;
00070 template<unsigned DIM> class VtkNonlinearElasticitySolutionWriter;
00071
00086 template<unsigned DIM>
00087 class AbstractContinuumMechanicsSolver
00088 {
00089 friend class StressRecoveror<DIM>;
00090 friend class VtkNonlinearElasticitySolutionWriter<DIM>;
00091
00092 protected:
00097 AbstractTetrahedralMesh<DIM, DIM>& mrQuadMesh;
00098
00100 ContinuumMechanicsProblemDefinition<DIM>& mrProblemDefinition;
00101
00103 bool mWriteOutput;
00104
00106 std::string mOutputDirectory;
00107
00109 OutputFileHandler* mpOutputFileHandler;
00110
00114 std::vector<c_vector<double,DIM> > mSpatialSolution;
00115
00117 std::vector<double> mPressureSolution;
00118
00119
00126 std::vector<double> mCurrentSolution;
00127
00129 GaussianQuadratureRule<DIM>* mpQuadratureRule;
00130
00132 GaussianQuadratureRule<DIM-1>* mpBoundaryQuadratureRule;
00133
00138 CompressibilityType mCompressibilityType;
00139
00148 unsigned mProblemDimension;
00149
00153 unsigned mNumDofs;
00154
00160 bool mVerbose;
00161
00182 Vec mResidualVector;
00183
00187 Vec mLinearSystemRhsVector;
00188
00192 Mat mSystemLhsMatrix;
00193
00197 Vec mDirichletBoundaryConditionsVector;
00198
00202 Mat mPreconditionMatrix;
00203
00207 void AllocateMatrixMemory();
00208
00209
00235 void ApplyDirichletBoundaryConditions(ApplyDirichletBcsType type, bool symmetricProblem);
00236
00269 void AddIdentityBlockForDummyPressureVariables(ApplyDirichletBcsType type);
00270
00271
00288 void RemovePressureDummyValuesThroughLinearInterpolation();
00289
00290 public:
00298 AbstractContinuumMechanicsSolver(AbstractTetrahedralMesh<DIM, DIM>& rQuadMesh,
00299 ContinuumMechanicsProblemDefinition<DIM>& rProblemDefinition,
00300 std::string outputDirectory,
00301 CompressibilityType compressibilityType);
00302
00304 virtual ~AbstractContinuumMechanicsSolver();
00305
00317 void WriteCurrentSpatialSolution(std::string fileName, std::string fileExtension, int counterToAppend=-1);
00318
00331 void WriteCurrentPressureSolution( int counterToAppend=-1);
00332
00338 void SetWriteOutput(bool writeOutput=true);
00339
00347 void CreateVtkOutput(std::string spatialSolutionName="Spatial solution");
00348
00353 std::vector<double>& rGetCurrentSolution()
00354 {
00355 return mCurrentSolution;
00356 }
00357
00362 virtual std::vector<c_vector<double,DIM> >& rGetSpatialSolution()=0;
00363
00372 std::vector<double>& rGetPressures();
00373
00374 };
00375
00376
00377 template<unsigned DIM>
00378 AbstractContinuumMechanicsSolver<DIM>::AbstractContinuumMechanicsSolver(AbstractTetrahedralMesh<DIM, DIM>& rQuadMesh,
00379 ContinuumMechanicsProblemDefinition<DIM>& rProblemDefinition,
00380 std::string outputDirectory,
00381 CompressibilityType compressibilityType)
00382 : mrQuadMesh(rQuadMesh),
00383 mrProblemDefinition(rProblemDefinition),
00384 mOutputDirectory(outputDirectory),
00385 mpOutputFileHandler(NULL),
00386 mpQuadratureRule(NULL),
00387 mpBoundaryQuadratureRule(NULL),
00388 mCompressibilityType(compressibilityType),
00389 mResidualVector(NULL),
00390 mSystemLhsMatrix(NULL),
00391 mPreconditionMatrix(NULL)
00392 {
00393 assert(DIM==2 || DIM==3);
00394
00395
00396 QuadraticMesh<DIM>* p_quad_mesh = dynamic_cast<QuadraticMesh<DIM>* >(&rQuadMesh);
00397 DistributedQuadraticMesh<DIM>* p_distributed_quad_mesh = dynamic_cast<DistributedQuadraticMesh<DIM>* >(&rQuadMesh);
00398
00399 if(p_quad_mesh == NULL && p_distributed_quad_mesh == NULL)
00400 {
00401 EXCEPTION("Continuum mechanics solvers require a quadratic mesh");
00402 }
00403
00404
00405 mVerbose = (mrProblemDefinition.GetVerboseDuringSolve() ||
00406 CommandLineArguments::Instance()->OptionExists("-mech_verbose") ||
00407 CommandLineArguments::Instance()->OptionExists("-mech_very_verbose") );
00408
00409 mWriteOutput = (mOutputDirectory != "");
00410 if (mWriteOutput)
00411 {
00412 mpOutputFileHandler = new OutputFileHandler(mOutputDirectory);
00413 }
00414
00415
00416 mProblemDimension = mCompressibilityType==COMPRESSIBLE ? DIM : DIM+1;
00417 mNumDofs = mProblemDimension*mrQuadMesh.GetNumNodes();
00418
00419 AllocateMatrixMemory();
00420
00421
00422
00423 mpQuadratureRule = new GaussianQuadratureRule<DIM>(3);
00424
00425
00426 mpBoundaryQuadratureRule = new GaussianQuadratureRule<DIM-1>(3);
00427
00428 mCurrentSolution.resize(mNumDofs, 0.0);
00429 }
00430
00431
00432 template<unsigned DIM>
00433 AbstractContinuumMechanicsSolver<DIM>::~AbstractContinuumMechanicsSolver()
00434 {
00435 if (mpOutputFileHandler)
00436 {
00437 delete mpOutputFileHandler;
00438 }
00439
00440 if (mpQuadratureRule)
00441 {
00442 delete mpQuadratureRule;
00443 delete mpBoundaryQuadratureRule;
00444 }
00445
00446 if (mResidualVector)
00447 {
00448 PetscTools::Destroy(mResidualVector);
00449 PetscTools::Destroy(mLinearSystemRhsVector);
00450 PetscTools::Destroy(mSystemLhsMatrix);
00451 PetscTools::Destroy(mPreconditionMatrix);
00452 }
00453
00454 if (mDirichletBoundaryConditionsVector)
00455 {
00456 PetscTools::Destroy(mDirichletBoundaryConditionsVector);
00457 }
00458 }
00459
00460 template<unsigned DIM>
00461 void AbstractContinuumMechanicsSolver<DIM>::WriteCurrentSpatialSolution(std::string fileName,
00462 std::string fileExtension,
00463 int counterToAppend)
00464 {
00465
00466 if (!mWriteOutput)
00467 {
00468 return;
00469 }
00470
00471 if (PetscTools::AmMaster())
00472 {
00473 std::stringstream file_name;
00474
00475 file_name << fileName;
00476 if (counterToAppend >= 0)
00477 {
00478 file_name << "_" << counterToAppend;
00479 }
00480 file_name << "." << fileExtension;
00481
00482 out_stream p_file = mpOutputFileHandler->OpenOutputFile(file_name.str());
00483
00484 std::vector<c_vector<double,DIM> >& r_spatial_solution = rGetSpatialSolution();
00485 for (unsigned i=0; i<r_spatial_solution.size(); i++)
00486 {
00487
00488
00489
00490
00491 for (unsigned j=0; j<DIM; j++)
00492 {
00493 *p_file << r_spatial_solution[i](j) << " ";
00494 }
00495 *p_file << "\n";
00496 }
00497 p_file->close();
00498 }
00499 PetscTools::Barrier("WriteSpatial");
00500 }
00501
00502 template<unsigned DIM>
00503 void AbstractContinuumMechanicsSolver<DIM>::WriteCurrentPressureSolution(int counterToAppend)
00504 {
00505
00506 if (!mWriteOutput)
00507 {
00508 return;
00509 }
00510
00511 if (PetscTools::AmMaster())
00512 {
00513 std::stringstream file_name;
00514
00515 file_name << "pressure";
00516 if (counterToAppend >= 0)
00517 {
00518 file_name << "_" << counterToAppend;
00519 }
00520 file_name << ".txt";
00521
00522 out_stream p_file = mpOutputFileHandler->OpenOutputFile(file_name.str());
00523
00524 std::vector<double>& r_pressure = rGetPressures();
00525 for (unsigned i=0; i<r_pressure.size(); i++)
00526 {
00527 for (unsigned j=0; j<DIM; j++)
00528 {
00529 *p_file << mrQuadMesh.GetNode(i)->rGetLocation()[j] << " ";
00530 }
00531
00532 *p_file << r_pressure[i] << "\n";
00533 }
00534 p_file->close();
00535 }
00536 PetscTools::Barrier("WritePressure");
00537 }
00538
00539
00540 template<unsigned DIM>
00541 void AbstractContinuumMechanicsSolver<DIM>::SetWriteOutput(bool writeOutput)
00542 {
00543 if (writeOutput && (mOutputDirectory==""))
00544 {
00545 EXCEPTION("Can't write output if no output directory was given in constructor");
00546 }
00547 mWriteOutput = writeOutput;
00548 }
00549
00550
00551 template<unsigned DIM>
00552 void AbstractContinuumMechanicsSolver<DIM>::CreateVtkOutput(std::string spatialSolutionName)
00553 {
00554 if (this->mOutputDirectory=="")
00555 {
00556 EXCEPTION("No output directory was given so no output was written, cannot convert to VTK format");
00557 }
00558 #ifdef CHASTE_VTK
00559 VtkMeshWriter<DIM, DIM> mesh_writer(this->mOutputDirectory + "/vtk", "solution", true);
00560
00561 mesh_writer.AddPointData(spatialSolutionName, this->rGetSpatialSolution());
00562
00563 if (mCompressibilityType==INCOMPRESSIBLE)
00564 {
00565 mesh_writer.AddPointData("Pressure", rGetPressures());
00566 }
00567
00568
00569 std::vector<double> element_attribute;
00570 for(typename QuadraticMesh<DIM>::ElementIterator iter = this->mrQuadMesh.GetElementIteratorBegin();
00571 iter != this->mrQuadMesh.GetElementIteratorEnd();
00572 ++iter)
00573 {
00574 element_attribute.push_back(iter->GetAttribute());
00575 }
00576 mesh_writer.AddCellData("Attribute", element_attribute);
00577
00578 mesh_writer.WriteFilesUsingMesh(this->mrQuadMesh);
00579 #endif
00580 }
00581
00582
00583 template<unsigned DIM>
00584 std::vector<double>& AbstractContinuumMechanicsSolver<DIM>::rGetPressures()
00585 {
00586 assert(mProblemDimension==DIM+1);
00587
00588 mPressureSolution.clear();
00589 mPressureSolution.resize(mrQuadMesh.GetNumNodes());
00590
00591 for (unsigned i=0; i<mrQuadMesh.GetNumNodes(); i++)
00592 {
00593 mPressureSolution[i] = mCurrentSolution[mProblemDimension*i + DIM];
00594 }
00595 return mPressureSolution;
00596 }
00597
00598
00599 template<unsigned DIM>
00600 void AbstractContinuumMechanicsSolver<DIM>::RemovePressureDummyValuesThroughLinearInterpolation()
00601 {
00602 assert(mProblemDimension==DIM+1);
00603
00604
00605 unsigned internal_nodes_2d[3] = {3,4,5};
00606 unsigned neighbouring_vertices_2d[3][2] = { {1,2}, {2,0}, {0,1} };
00607
00608
00609 unsigned internal_nodes_3d[6] = {4,5,6,7,8,9};
00610 unsigned neighbouring_vertices_3d[6][2] = { {0,1}, {1,2}, {0,2}, {0,3}, {1,3}, {2,3} };
00611
00612 unsigned num_internal_nodes_per_element = DIM==2 ? 3 : 6;
00613
00614
00615 for (typename AbstractTetrahedralMesh<DIM,DIM>::ElementIterator iter = mrQuadMesh.GetElementIteratorBegin();
00616 iter != mrQuadMesh.GetElementIteratorEnd();
00617 ++iter)
00618 {
00619 for(unsigned i=0; i<num_internal_nodes_per_element; i++)
00620 {
00621 unsigned global_index;
00622 double left_val;
00623 double right_val;
00624
00625 if(DIM==2)
00626 {
00627 global_index = iter->GetNodeGlobalIndex( internal_nodes_2d[i] );
00628 unsigned vertex_0_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_2d[i][0] );
00629 unsigned vertex_1_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_2d[i][1] );
00630 left_val = mCurrentSolution[mProblemDimension*vertex_0_global_index + DIM];
00631 right_val = mCurrentSolution[mProblemDimension*vertex_1_global_index + DIM];
00632 }
00633 else
00634 {
00635 global_index = iter->GetNodeGlobalIndex( internal_nodes_3d[i] );
00636 unsigned vertex_0_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_3d[i][0] );
00637 unsigned vertex_1_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_3d[i][1] );
00638 left_val = mCurrentSolution[mProblemDimension*vertex_0_global_index + DIM];
00639 right_val = mCurrentSolution[mProblemDimension*vertex_1_global_index + DIM];
00640 }
00641
00642
00643 mCurrentSolution[mProblemDimension*global_index + DIM] = 0.5 * (left_val + right_val);
00644 }
00645 }
00646 }
00647
00648
00649
00650
00651
00652
00653
00654
00655
00656
00657
00658
00659
00660
00661
00662
00663
00664
00665
00666
00667
00668
00669
00670
00671
00672
00673
00674
00675
00676
00677
00678 template<unsigned DIM>
00679 void AbstractContinuumMechanicsSolver<DIM>::ApplyDirichletBoundaryConditions(ApplyDirichletBcsType type, bool symmetricProblem)
00680 {
00681 std::vector<unsigned> rows;
00682 std::vector<double> values;
00683
00684
00685 bool applySymmetrically = (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY) && symmetricProblem;
00686
00687 if (applySymmetrically)
00688 {
00689 if(mDirichletBoundaryConditionsVector==NULL)
00690 {
00691 VecDuplicate(mResidualVector, &mDirichletBoundaryConditionsVector);
00692 }
00693
00694 PetscVecTools::Zero(mDirichletBoundaryConditionsVector);
00695 PetscMatTools::Finalise(mSystemLhsMatrix);
00696 }
00697
00699
00701
00702 for (unsigned i=0; i<mrProblemDefinition.rGetDirichletNodes().size(); i++)
00703 {
00704 unsigned node_index = mrProblemDefinition.rGetDirichletNodes()[i];
00705
00706 for (unsigned j=0; j<DIM; j++)
00707 {
00708 double dirichlet_val = mrProblemDefinition.rGetDirichletNodeValues()[i](j);
00709
00710 if(dirichlet_val != ContinuumMechanicsProblemDefinition<DIM>::FREE)
00711 {
00712 double val;
00713 unsigned dof_index = mProblemDimension*node_index+j;
00714
00715 if(type == LINEAR_PROBLEM)
00716 {
00717 val = dirichlet_val;
00718 }
00719 else
00720 {
00721
00722
00723
00724
00725
00726 val = mCurrentSolution[dof_index] - dirichlet_val;
00727 }
00728 rows.push_back(dof_index);
00729 values.push_back(val);
00730 }
00731 }
00732 }
00733
00735
00737
00738 if (applySymmetrically)
00739 {
00740
00741 for (unsigned i=0; i<rows.size(); i++)
00742 {
00743 unsigned col = rows[i];
00744 double minus_value = -values[i];
00745
00746
00747
00748
00749
00750
00751 Vec matrix_col = PetscMatTools::GetMatrixRowDistributed(mSystemLhsMatrix,col);
00752
00753
00754 PetscVecTools::SetElement(matrix_col, col, 0.0);
00755
00756
00757
00758
00759
00760 PetscVecTools::AddScaledVector(mDirichletBoundaryConditionsVector, matrix_col, minus_value);
00761 PetscTools::Destroy(matrix_col);
00762 }
00763 }
00764
00765 if (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY)
00766 {
00767
00768
00769
00770 if (applySymmetrically)
00771 {
00772 PetscMatTools::ZeroRowsAndColumnsWithValueOnDiagonal(mSystemLhsMatrix, rows, 1.0);
00773 PetscMatTools::ZeroRowsAndColumnsWithValueOnDiagonal(mPreconditionMatrix, rows, 1.0);
00774
00775
00776 PetscVecTools::AddScaledVector(mLinearSystemRhsVector, mDirichletBoundaryConditionsVector, 1.0);
00777 }
00778 else
00779 {
00780 PetscMatTools::ZeroRowsWithValueOnDiagonal(mSystemLhsMatrix, rows, 1.0);
00781 PetscMatTools::ZeroRowsWithValueOnDiagonal(mPreconditionMatrix, rows, 1.0);
00782 }
00783 }
00784
00785 if (type!=LINEAR_PROBLEM)
00786 {
00787 for (unsigned i=0; i<rows.size(); i++)
00788 {
00789 PetscVecTools::SetElement(mResidualVector, rows[i], values[i]);
00790 }
00791 }
00792
00793 if (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY)
00794 {
00795 for (unsigned i=0; i<rows.size(); i++)
00796 {
00797 PetscVecTools::SetElement(mLinearSystemRhsVector, rows[i], values[i]);
00798 }
00799 }
00800 }
00801
00802
00803 template<unsigned DIM>
00804 void AbstractContinuumMechanicsSolver<DIM>::AddIdentityBlockForDummyPressureVariables(ApplyDirichletBcsType type)
00805 {
00806 assert(mCompressibilityType==INCOMPRESSIBLE);
00807
00808 int lo, hi;
00809 VecGetOwnershipRange(mResidualVector, &lo, &hi);
00810
00811 for(unsigned i=0; i<mrQuadMesh.GetNumNodes(); i++)
00812 {
00813 if(mrQuadMesh.GetNode(i)->IsInternal())
00814 {
00815 unsigned row = (DIM+1)*i + DIM;
00816 if(lo <= (int)row && (int)row < hi)
00817 {
00818 if(type!=LINEAR_PROBLEM)
00819 {
00820 PetscVecTools::SetElement(mResidualVector, row, mCurrentSolution[row]-0.0);
00821 }
00822 if(type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY)
00823 {
00824 double rhs_vector_val = type==LINEAR_PROBLEM ? 0.0 : mCurrentSolution[row]-0.0;
00825 PetscVecTools::SetElement(mLinearSystemRhsVector, row, rhs_vector_val);
00826
00827 PetscMatTools::SetElement(mSystemLhsMatrix, row, row, 1.0);
00828 PetscMatTools::SetElement(mPreconditionMatrix, row, row, 1.0);
00829 }
00830 }
00831 }
00832 }
00833 }
00834
00835
00836
00837 template<unsigned DIM>
00838 void AbstractContinuumMechanicsSolver<DIM>::AllocateMatrixMemory()
00839 {
00840 Vec template_vec = mrQuadMesh.GetDistributedVectorFactory()->CreateVec(mProblemDimension);
00841
00843
00845 VecDuplicate(template_vec, &mResidualVector);
00846 VecDuplicate(mResidualVector, &mLinearSystemRhsVector);
00847
00848
00849 mDirichletBoundaryConditionsVector = NULL;
00850 PetscTools::Destroy(template_vec);
00851
00853
00855
00856 int lo, hi;
00857 VecGetOwnershipRange(mResidualVector, &lo, &hi);
00858 PetscInt local_size = hi - lo;
00859
00860
00861 if (DIM==2)
00862 {
00863
00864 unsigned num_non_zeros = std::min(75u, mNumDofs);
00865
00866 PetscTools::SetupMat(mSystemLhsMatrix, mNumDofs, mNumDofs, num_non_zeros, local_size, local_size);
00867 PetscTools::SetupMat(mPreconditionMatrix, mNumDofs, mNumDofs, num_non_zeros, local_size, local_size);
00868 }
00869 else
00870 {
00871 assert(DIM==3);
00872
00873
00874
00875
00876 int* num_non_zeros_each_row = new int[mNumDofs];
00877 for (unsigned i=0; i<mNumDofs; i++)
00878 {
00879 num_non_zeros_each_row[i] = 0;
00880 }
00881
00882 for (typename AbstractMesh<DIM,DIM>::NodeIterator iter = mrQuadMesh.GetNodeIteratorBegin();
00883 iter != mrQuadMesh.GetNodeIteratorEnd();
00884 ++iter)
00885 {
00886
00887
00888
00889 unsigned num_non_zeros_upper_bound = 4 + 30*iter->GetNumContainingElements();
00890
00891 num_non_zeros_upper_bound = std::min(num_non_zeros_upper_bound, mNumDofs);
00892
00893 unsigned i = iter->GetIndex();
00894
00895 num_non_zeros_each_row[mProblemDimension*i + 0] = num_non_zeros_upper_bound;
00896 num_non_zeros_each_row[mProblemDimension*i + 1] = num_non_zeros_upper_bound;
00897 num_non_zeros_each_row[mProblemDimension*i + 2] = num_non_zeros_upper_bound;
00898
00899 if (mCompressibilityType==INCOMPRESSIBLE)
00900 {
00901 if(!iter->IsInternal())
00902 {
00903 num_non_zeros_each_row[mProblemDimension*i + 3] = num_non_zeros_upper_bound;
00904 }
00905 else
00906 {
00907 num_non_zeros_each_row[mProblemDimension*i + 3] = 1;
00908 }
00909 }
00910 }
00911
00912
00913
00914
00915
00916
00919
00920
00922
00923
00924
00925
00926 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00927 MatCreate(PETSC_COMM_WORLD,local_size,local_size,mNumDofs,mNumDofs,&mSystemLhsMatrix);
00928 MatCreate(PETSC_COMM_WORLD,local_size,local_size,mNumDofs,mNumDofs,&mPreconditionMatrix);
00929 #else //New API
00930 MatCreate(PETSC_COMM_WORLD,&mSystemLhsMatrix);
00931 MatCreate(PETSC_COMM_WORLD,&mPreconditionMatrix);
00932 MatSetSizes(mSystemLhsMatrix,local_size,local_size,mNumDofs,mNumDofs);
00933 MatSetSizes(mPreconditionMatrix,local_size,local_size,mNumDofs,mNumDofs);
00934 #endif
00935
00936 if (PetscTools::IsSequential())
00937 {
00938 MatSetType(mSystemLhsMatrix, MATSEQAIJ);
00939 MatSetType(mPreconditionMatrix, MATSEQAIJ);
00940 MatSeqAIJSetPreallocation(mSystemLhsMatrix, PETSC_DEFAULT, num_non_zeros_each_row);
00941 MatSeqAIJSetPreallocation(mPreconditionMatrix, PETSC_DEFAULT, num_non_zeros_each_row);
00942 }
00943 else
00944 {
00945 int* num_non_zeros_each_row_in_diag = new int[local_size];
00946 int* num_non_zeros_each_row_off_diag = new int[local_size];
00947 for (unsigned i=0; i<unsigned(local_size); i++)
00948 {
00949 num_non_zeros_each_row_in_diag[i] = num_non_zeros_each_row[lo+i];
00950 num_non_zeros_each_row_off_diag[i] = num_non_zeros_each_row[lo+i];
00951
00952 if(num_non_zeros_each_row_in_diag[i] > local_size)
00953 {
00954 num_non_zeros_each_row_in_diag[i] = local_size;
00955 }
00956 }
00957
00958 MatSetType(mSystemLhsMatrix, MATMPIAIJ);
00959 MatSetType(mPreconditionMatrix, MATMPIAIJ);
00960 MatMPIAIJSetPreallocation(mSystemLhsMatrix, PETSC_DEFAULT, num_non_zeros_each_row_in_diag, PETSC_DEFAULT, num_non_zeros_each_row_off_diag);
00961 MatMPIAIJSetPreallocation(mPreconditionMatrix, PETSC_DEFAULT, num_non_zeros_each_row_in_diag, PETSC_DEFAULT, num_non_zeros_each_row_off_diag);
00962 }
00963
00964 MatSetFromOptions(mSystemLhsMatrix);
00965 MatSetFromOptions(mPreconditionMatrix);
00966 #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
00967 MatSetOption(mSystemLhsMatrix, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE);
00968 MatSetOption(mPreconditionMatrix, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE);
00969 #else
00970 MatSetOption(mSystemLhsMatrix, MAT_IGNORE_OFF_PROC_ENTRIES);
00971 MatSetOption(mPreconditionMatrix, MAT_IGNORE_OFF_PROC_ENTRIES);
00972 #endif
00973
00974
00975
00976
00977
00978
00979
00980
00981 delete [] num_non_zeros_each_row;
00982 }
00983 }
00984 #endif // ABSTRACTCONTINUUMMECHANICSSOLVER_HPP_