svFSIplus
|
wrap Trilinos solver functions More...
#include <stdio.h>
#include <vector>
#include <iostream>
#include "mpi.h"
#include <time.h>
#include <numeric>
#include "Epetra_MpiComm.h"
#include "Epetra_Map.h"
#include "Epetra_FEVbrMatrix.h"
#include "Epetra_FEVector.h"
#include "Epetra_FECrsGraph.h"
#include "Epetra_FECrsMatrix.h"
#include "Epetra_Import.h"
#include "AztecOO.h"
#include "AztecOO_StatusTestResNorm.h"
#include "ml_include.h"
#include "ml_MultiLevelPreconditioner.h"
#include "Ifpack.h"
#include "Ifpack_ConfigDefs.h"
#include "Ifpack_AdditiveSchwarz.h"
#include "Ifpack_ILUT.h"
#include "Ifpack_IC.h"
#include "Ifpack_ICT.h"
Go to the source code of this file.
Classes | |
struct | Trilinos |
Initialize all Epetra types we need separate from Fortran. More... | |
class | TrilinosMatVec |
This class implements the pure virtual class Epetra_Operator for the AztecOO iterative solve which only uses the Apply() method to compute the matrix vector product. More... | |
Functions | |
void | trilinos_lhs_create_ (int &numGlobalNodes, int &numLocalNodes, int &numGhostAndLocalNodes, int &nnz, const int *ltgSorted, const int *ltgUnsorted, const int *rowPtr, const int *colInd, int &dof, int &cpp_index, int &proc_id) |
Give function definitions which will be called through fortran. More... | |
void | trilinos_bc_create_ (const double *v, bool &isCoupledBC) |
void | trilinos_doassem_ (int &numNodesPerElement, const int *eqN, const double *lK, double *lR) |
void | trilinos_global_solve_ (const double *Val, const double *RHS, double *x, const double *dirW, double &resNorm, double &initNorm, int &numIters, double &solverTime, double &dB, bool &converged, int &lsType, double &relTol, int &maxIters, int &kspace, int &precondType) |
void | trilinos_solve_ (double *x, const double *dirW, double &resNorm, double &initNorm, int &numIters, double &solverTime, double &dB, bool &converged, int &lsType, double &relTol, int &maxIters, int &kspace, int &precondType, bool &isFassem) |
void | trilinos_lhs_free_ () |
void | setPreconditioner (int precondType, AztecOO &Solver) |
void | setMLPrec (AztecOO &Solver) |
void | setIFPACKPrec (AztecOO &Solver) |
void | checkDiagonalIsZero () |
void | constructJacobiScaling (const double *dirW, Epetra_Vector &diagonal) |
void | printMatrixToFile () |
void | printRHSToFile () |
void | printSolutionToFile () |
wrap Trilinos solver functions
Copyright (c) Stanford University, The Regents of the University of California, and others.
All Rights Reserved.
See Copyright-SimVascular.txt for additional details.
Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the "Software"), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions:
The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
void checkDiagonalIsZero | ( | ) |
This routine is to be used with preconditioners such as ILUT which require 1s on the diagonal
void constructJacobiScaling | ( | const double * | dirW, |
Epetra_Vector & | diagonal | ||
) |
To be called within solve-output diagonal from here
dirW | pass in array with Dirichlet boundary face nodes marked \paramdiagonal diagonal scaling vector need to output to multiply solution by |
void printMatrixToFile | ( | ) |
for debugging purposes here are routines to print the matrix and RHS vector
void setIFPACKPrec | ( | AztecOO & | Solver | ) |
pass in IC, ICT pass in string for which to turn on right now set to IC
void setMLPrec | ( | AztecOO & | Solver | ) |
Tune parameters for htis and IFPACK Ref: https://trilinos.org/oldsite/packages/ml/mlguide5.pdf
void trilinos_bc_create_ | ( | const double * | v, |
bool & | isCoupledBC | ||
) |
v | coeff in the scalar product |
isCoupledBC | determines if coupled resistance BC is turned on |
v | coupled boundary vector |
isCoupledBC | determines if coupled resistance BC is turned on |
void trilinos_doassem_ | ( | int & | numNodesPerElement, |
const int * | eqN, | ||
const double * | lK, | ||
double * | lR | ||
) |
This function assembles the dense element block stiffness matrix and element block force vector into the global block stiffness K and global block force vector F. Happens on a single element and will be called repeatedly as we loop over all the elements.
numNodesPerElement | number of nodes (d) for element e |
eqN | converts local element node number to local proc equation number-size is numNodesPerElement |
lK | element local stiffness of size (dof*dof, numNodesPerElement, numNodesPerElement) |
lR | element local force vector of size (dof, numNodesPerElement) |
void trilinos_global_solve_ | ( | const double * | Val, |
const double * | RHS, | ||
double * | x, | ||
const double * | dirW, | ||
double & | resNorm, | ||
double & | initNorm, | ||
int & | numIters, | ||
double & | solverTime, | ||
double & | dB, | ||
bool & | converged, | ||
int & | lsType, | ||
double & | relTol, | ||
int & | maxIters, | ||
int & | kspace, | ||
int & | precondType | ||
) |
Tthis function creates the LHS structure topology and RHS vector based on the block map using the global stiffness matrix and global force vector assuming the assembly part has been done in svFSI
Val | nonzero values of assembled global stiffness |
RHS | global force vector |
x | solution vector |
dirW | gives information about the Dirichlet BC |
resNorm | norm of solution |
initNorm | initial norm of RHS |
numIters | linear solver number of iterations |
solverTime | time in linear solver |
dB | log ratio for output |
converged | neither or not converged in the max number of iterations |
lsType | defines type of linear solver to use |
relTol | default relative tolerance or as set in input file |
maxIters | default max number of iterations for gmres per restart |
kspace | specific for gmres dim of the stored Krylov space vectors |
precondType | defines type of preconditioner to use |
void trilinos_lhs_create_ | ( | int & | numGlobalNodes, |
int & | numLocalNodes, | ||
int & | numGhostAndLocalNodes, | ||
int & | nnz, | ||
const int * | ltgSorted, | ||
const int * | ltgUnsorted, | ||
const int * | rowPtr, | ||
const int * | colInd, | ||
int & | Dof, | ||
int & | cpp_index, | ||
int & | proc_id | ||
) |
Give function definitions which will be called through fortran.
make the graph global and optimize storage on it this function creates the LHS structure topology and RHS vector based on the block map
numGlobalNodes | total/global number of nodes each node can have dof for the spatial dim |
numLocalNodes | number of nodes owned by this proc in block coordiantes |
numGhostAndLocalNodes | number of nodes owned and shared by this processor includes ghost nodes |
nnz | number of nonzeros in LHS matrix of calling proc |
ltgSorted | integer pointer of size numLocalNodes returns global node number of local node a to give blockrow |
ltgUnsorted | unsorted/not-reordered version |
rowPtr | CSR row ptr of size numLocalNodes + 1 to block rows |
colInd | CSR column indices ptr (size nnz points) to block rows |
Dof | size of each block element to give dim of each block |
void trilinos_lhs_free_ | ( | ) |
free the memory of the global data structures
void trilinos_solve_ | ( | double * | x, |
const double * | dirW, | ||
double & | resNorm, | ||
double & | initNorm, | ||
int & | numIters, | ||
double & | solverTime, | ||
double & | dB, | ||
bool & | converged, | ||
int & | lsType, | ||
double & | relTol, | ||
int & | maxIters, | ||
int & | kspace, | ||
int & | precondType, | ||
bool & | isFassem | ||
) |
This function uses the established data structure to solve the block linear system and passes the solution vector back to fortran with the ghost nodes filled in
x | solution vector |
dirW | gives information about Dirichlet BC |
resNorm | norm of the final residual |
initNorm | norm of initial resiual x_init = 0 ||b|| |
numIters | number of iterations in the linear solver |
solverTime | time in the linear solver |
dB | log ratio for output |
converged | can pass in solver and preconditioner type too |
lsType | defines type of linear solver to use |
relTol | default relative tolerance or as set in input file |
maxIters | default max number of iterations for gmres per restart |
kspace | specific for gmres dim of the stored Krylov space vectors |
precondType | defines type of preconditioner to use |
isFassem | determines if F is already assembled at ghost nodes |