svFSIplus
petsc_impl.h
1 /* Copyright (c) Stanford University, The Regents of the University of California, and others.
2  *
3  * All Rights Reserved.
4  *
5  * See Copyright-SimVascular.txt for additional details.
6  *
7  * Permission is hereby granted, free of charge, to any person obtaining
8  * a copy of this software and associated documentation files (the
9  * "Software"), to deal in the Software without restriction, including
10  * without limitation the rights to use, copy, modify, merge, publish,
11  * distribute, sublicense, and/or sell copies of the Software, and to
12  * permit persons to whom the Software is furnished to do so, subject
13  * to the following conditions:
14  *
15  * The above copyright notice and this permission notice shall be included
16  * in all copies or substantial portions of the Software.
17  *
18  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS
19  * IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
20  * TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A
21  * PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER
22  * OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
23  * EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
24  * PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
25  * PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
26  * LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
27  * NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
28  * SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
29  */
30 
31 // This file contains PETSc-dependent data structures.
32 
33 #ifndef PETSC_INTERFACE_H
34 #define PETSC_INTERFACE_H
35 
36 #include <petscksp.h>
37 #include <petscao.h>
38 #include <unistd.h>
39 #include <stdbool.h>
40 
41 #include "consts.h"
42 
43 //--------
44 // LHSCtx
45 //---------
46 // PETSc lhs context.
47 //
48 typedef struct {
49  PetscBool created; /* Whether petsc lhs is created */
50 
51  PetscInt nNo; /* local number of vertices */
52  PetscInt mynNo; /* number of owned vertices */
53 
54  PetscInt *map; /* local to local mapping, map[O2] = O1 */
55  PetscInt *ltg; /* local to global in PETSc ordering */
56  PetscInt *ghostltg; /* local to global in PETSc ordering */
57  PetscInt *rowPtr; /* row pointer for adjacency info */
58  PetscInt *colPtr; /* column pointer for adjacency info */
59 } LHSCtx;
60 
61 //-------
62 // LSCtx
63 //-------
64 // PETSc linear solver context.
65 //
66 typedef struct {
67  PetscBool created; /* Whether mat and vec is created */
68  const char *pre; /* option prefix for different equations */
69 
70  PetscInt lpPts; /* number of dofs with lumped parameter BC */
71  PetscInt *lpBC_l; /* O2 index for dofs with lumped parameter BC */
72  PetscInt *lpBC_g; /* PETSc index for dofs with lumped parameter BC */
73 
74  PetscInt DirPts; /* number of dofs with Dirichlet BC */
75  PetscInt *DirBC; /* PETSc index for dofs with Dirichlet BC */
76 
77  Vec b; /* rhs/solution vector of owned vertices */
78  Mat A; /* stiffness matrix */
79  KSP ksp; /* linear solver context */
80 
81  PetscBool rcs; /* whether rcs preconditioner is activated */
82  Vec Dr; /* diagonal matrix from row maxabs */
83  Vec Dc; /* diagonal matrix from col maxabs */
84 } LSCtx;
85 
86 void petsc_initialize(const PetscInt nNo, const PetscInt mynNo,
87  const PetscInt nnz, const PetscInt nEq, const PetscInt *svFSI_ltg,
88  const PetscInt *svFSI_map, const PetscInt *svFSI_rowPtr,
89  const PetscInt *svFSI_colPtr, char *inp);
90 
91 void petsc_create_linearsystem(const PetscInt dof, const PetscInt iEq, const PetscInt nEq,
92  const PetscReal *svFSI_DirBC, const PetscReal *svFSI_lpBC);
93 
94 void petsc_create_linearsolver(const consts::PreconditionerType lsType, const consts::PreconditionerType pcType,
95  const PetscInt kSpace, const PetscInt maxIter, const PetscReal relTol,
96  const PetscReal absTol, const consts::EquationType phys, const PetscInt dof,
97  const PetscInt iEq, const PetscInt nEq);
98 
99 void petsc_set_values(const PetscInt dof, const PetscInt iEq, const PetscReal *R,
100  const PetscReal *Val, const PetscReal *svFSI_DirBC, const PetscReal *svFSI_lpBC);
101 
102 void petsc_solve(PetscReal *resNorm, PetscReal *initNorm, PetscReal *dB,
103  PetscReal *execTime, bool *converged, PetscInt *numIter,
104  PetscReal *R, const PetscInt maxIter, const PetscInt dof,
105  const PetscInt iEq);
106 
107 void petsc_destroy_all(const PetscInt);
108 
109 PetscErrorCode petsc_create_lhs(const PetscInt, const PetscInt, const PetscInt,
110  const PetscInt *, const PetscInt *,
111  const PetscInt *, const PetscInt *);
112 
113 PetscErrorCode petsc_create_bc(const PetscInt, const PetscInt, const PetscReal *,
114  const PetscReal *);
115 
116 PetscErrorCode petsc_create_vecmat(const PetscInt, const PetscInt, const PetscInt);
117 
118 PetscErrorCode petsc_set_vec(const PetscInt, const PetscInt, const PetscReal *);
119 
120 PetscErrorCode petsc_set_mat(const PetscInt, const PetscInt, const PetscReal *);
121 
122 PetscErrorCode petsc_set_bc(const PetscInt, const PetscReal *, const PetscReal *);
123 
124 PetscErrorCode petsc_set_pcfieldsplit(const PetscInt, const PetscInt);
125 
126 PetscErrorCode petsc_pc_rcs(const PetscInt, const PetscInt);
127 
128 
129 PetscErrorCode petsc_debug_save_vec(const char *, Vec);
130 PetscErrorCode petsc_debug_save_mat(const char *, Mat);
131 
132 #endif
Linear system of equations solver type.
Definition: ComMod.h:640
Definition: petsc_impl.h:48
Definition: petsc_impl.h:66