Go to the documentation of this file.
14 #include "runTimeSelectionTables.H"
15 #include "fvOptions.H"
21 #include "syncTools.H"
109 const label connectedLevelLocal,
110 const wordList connectedStatesLocal,
111 const List<List<word>> connectedStateInterProc,
112 const label addFace);
117 const label idx)
const;
127 Mat* stateBoundaryCon,
128 Mat* stateBoundaryConTmp);
135 const word stateName,
136 const PetscScalar val);
143 const word stateName,
144 const PetscScalar val);
151 const word stateName,
152 const PetscScalar val);
160 const List<List<word>> connectedStates,
161 const label addFaces);
188 const word stateName,
190 const label comp = -1);
201 (
const word modelType,
216 const word modelType,
225 static autoPtr<DAJacCon>
New(
226 const word modelType,
240 virtual void clear() = 0;
249 virtual void setupJacCon(
const dictionary& options) = 0;
260 const label transposed)
const;
280 Vec coloredColumn)
const;
284 scalarList objFuncFaceValues,
285 scalarList objFuncCellValues,
286 Vec objFuncVec)
const;
const word modelType_
the name of the jacCon matrix
DAColoring daColoring_
DAColoring object.
void setPIVVec(Vec iSPIV, const label idxI)
function used to add connectivity for pressureInletVelocity
void checkSpecialBCs()
check if there is special boundary conditions that need special treatment in jacCon_
void calcColoredColumns(const label colorI, Vec coloredColumn) const
calculate the colored column vector
const fvMesh & mesh_
fvMesh
Mat jacCon_
Jacobian connectivity mat.
virtual void clear()=0
clear members in parent and child objects
const DAOption & daOption_
DAOption object.
void clearDAJacConMembers()
clear members in DAJacCon
label addPhi4PIV(const word stateName, const label idxI, const label comp=-1)
add connectivity phi for pressureInletVelocity
virtual void setObjFuncVec(scalarList objFuncFaceValues, scalarList objFuncCellValues, Vec objFuncVec) const
assign values for the objective function vector based on the face and cell value lists
DAOption daOption(mesh, pyOptions_)
void setupStateBoundaryConID(Mat *stateBoundaryConID)
calculate DAJacCon::stateBoundaryConID_
Vec getJacConColor() const
return DAJacCon::jacConColors_
void combineStateBndCon(Mat *stateBoundaryCon, Mat *stateBoundaryConTmp)
combine stateBoundaryCon and stateBoundaryConTmp, this is to ensure including all connected states fo...
label getLocalCoupledBFaceIndex(const label localFaceI) const
given a local face index, return the local index of the coupled boundary face
TypeName("DAJacCon")
Runtime type information.
label getNJacConColors() const
get the number of JacCon colors
const DAModel & daModel_
DAModel object.
void initializeStateBoundaryCon()
initialize state boundary connection
void setupJacobianConnections(Mat conMat, Mat connections, const label idxI)
assign values in connections to a specific row idxI in conMat
void calcNeiBFaceGlobalCompact(labelList &neiBFaceGlobalCompact)
calculate DAJacCon::neiBFaceGlobalCompact_
void addConMatCell(Mat conMat, const label gRow, const label cellI, const word stateName, const PetscScalar val)
add val to the gRow row in conMat, the column indice are the state (given by stateName) for a given c...
void addConMatNeighbourCells(Mat conMat, const label gRow, const label cellI, const word stateName, const PetscScalar val)
add val to gRow row in conMat, column indice are the neighbouring states (given by stateName) for a g...
virtual void setupJacCon(const dictionary &options)=0
assign 1 to all non-zero elements for the Jacobian connecitivyt matrix
virtual void preallocatedRdW(Mat dRMat, const label transposed) const
preallocate dRdW matrix using the preallocVec
Vec isPIVBCState_
a vector to show whether a state is connected to a pressureInletVelocity boundary face (3 level max)
virtual void initializeJacCon(const dictionary &options)=0
initialize the state Jacobian connectivity matrix
Mat stateBoundaryCon_
matrix to store boundary connectivity levels for state Jacobians
DAField daField_
DAField object.
const DAIndex & daIndex_
DAIndex object.
label coloringExists(const word postFix="") const
whether the coloring file exists
void setupStateBoundaryCon(Mat *stateBoundaryCon)
calculate DAJacCon::stateBoundaryCon_
label nJacConColors_
number of jacCon colors
labelList neiBFaceGlobalCompact_
neibough face global index for a given local boundary face
DAModel daModel(mesh, daOption)
void createConnectionMat(Mat *connectedStates)
allocate connectedState matrix
virtual void setupJacConPreallocation(const dictionary &options)
calculate the preallocation vector for initializing the JacCon mat, if necessary
void setConnections(Mat conMat, const label idx) const
add value 1 for the colume idx to conMat
void calcJacConColoring(const word postFix="")
compute graph coloring for Jacobian connectivity matrix
Mat stateBoundaryConID_
matrix to store boundary connectivity ID for state Jacobians
DAIndex daIndex(mesh, daOption, daModel)
Vec jacConColors_
jacCon matrix colors
void addConMatCellFaces(Mat conMat, const label gRow, const label cellI, const word stateName, const PetscScalar val)
add val to gRow row in conMat, column indice are the faces (given by stateName) for a given cellI
void addBoundaryFaceConnections(Mat conMat, const label gRow, const label cellI, const labelList v, const List< List< word >> connectedStates, const label addFaces)
add the column index of the (iner-proc) connected states and faces to conMat, given a local face inde...
void readJacConColoring(const word postFix="")
read colors for JacCon
void addStateConnections(Mat connections, const label cellI, const label connectedLevelLocal, const wordList connectedStatesLocal, const List< List< word >> connectedStateInterProc, const label addFace)
a high-level function to add connected state column indices to the connectivity matrix
HashTable< wordList > stateInfo_
the regState_ list from DAStateInfo object
declareRunTimeSelectionTable(autoPtr, DAJacCon, dictionary,(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex),(modelType, mesh, daOption, daModel, daIndex))
static autoPtr< DAJacCon > New(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex)