14 #include "../config/config.hpp"
20 #include "../fem/fem.hpp"
23 #if defined(PETSC_HAVE_HYPRE)
24 #include "petscmathypre.h"
28 #if PETSC_VERSION_LT(3,12,0)
29 #define MatComputeOperator(A,B,C) MatComputeExplicitOperator(A,C)
30 #define MatComputeOperatorTranspose(A,B,C) MatComputeExplicitOperatorTranspose(A,C)
43 #define PCHKERRQ(obj,err) do { \
46 PetscError(PetscObjectComm((PetscObject)(obj)),__LINE__,_MFEM_FUNC_NAME, \
47 __FILE__,(err),PETSC_ERROR_REPEAT,NULL); \
48 MFEM_ABORT("Error in PETSc. See stacktrace above."); \
51 #define CCHKERRQ(comm,err) do { \
54 PetscError(comm,__LINE__,_MFEM_FUNC_NAME, \
55 __FILE__,(err),PETSC_ERROR_REPEAT,NULL); \
56 MFEM_ABORT("Error in PETSc. See stacktrace above."); \
72 static PetscErrorCode __mfem_snes_jacobian(
SNES,
Vec,
Mat,
Mat,
void*);
73 static PetscErrorCode __mfem_snes_function(
SNES,
Vec,
Vec,
void*);
76 static PetscErrorCode __mfem_snes_postcheck(SNESLineSearch,
Vec,
Vec,
Vec,
77 PetscBool*,PetscBool*,
void*);
79 static PetscErrorCode __mfem_pc_shell_apply(
PC,
Vec,
Vec);
80 static PetscErrorCode __mfem_pc_shell_apply_transpose(
PC,
Vec,
Vec);
81 static PetscErrorCode __mfem_pc_shell_setup(
PC);
82 static PetscErrorCode __mfem_pc_shell_destroy(
PC);
83 static PetscErrorCode __mfem_pc_shell_view(
PC,PetscViewer);
84 static PetscErrorCode __mfem_mat_shell_apply(
Mat,
Vec,
Vec);
85 static PetscErrorCode __mfem_mat_shell_apply_transpose(
Mat,
Vec,
Vec);
86 static PetscErrorCode __mfem_mat_shell_destroy(
Mat);
87 static PetscErrorCode __mfem_mat_shell_copy(
Mat,
Mat,MatStructure);
88 static PetscErrorCode __mfem_array_container_destroy(
void*);
89 static PetscErrorCode __mfem_matarray_container_destroy(
void*);
90 static PetscErrorCode __mfem_monitor_ctx_destroy(
void**);
93 static PetscErrorCode Convert_Array_IS(MPI_Comm,
bool,
const mfem::Array<int>*,
98 static PetscErrorCode MakeShellPCWithFactory(
PC,
104 #if !defined(PETSC_HAVE_HYPRE)
105 static PetscErrorCode MatConvert_hypreParCSR_AIJ(hypre_ParCSRMatrix*,
Mat*);
106 static PetscErrorCode MatConvert_hypreParCSR_IS(hypre_ParCSRMatrix*,
Mat*);
115 unsigned long int numprec;
116 } __mfem_pc_shell_ctx;
127 void (*postcheck)(
mfem::Operator *op,
const mfem::Vector&, mfem::Vector&,
128 mfem::Vector&,
bool&,
bool&);
132 const mfem::Vector&,
const mfem::Vector&,
133 const mfem::Vector&,
const mfem::Vector&);
145 PetscObjectState cached_ijacstate;
146 PetscObjectState cached_rhsjacstate;
147 PetscObjectState cached_splits_xstate;
148 PetscObjectState cached_splits_xdotstate;
155 } __mfem_monitor_ctx;
158 static PetscErrorCode ierr;
178 ierr = PetscInitialize(argc,argv,rc_file,help);
179 MFEM_VERIFY(!ierr,
"Unable to initialize PETSc");
184 ierr = PetscFinalize();
185 MFEM_VERIFY(!ierr,
"Unable to finalize PETSc");
190 void PetscParVector::_SetDataAndSize_()
195 ierr = VecGetArrayRead(x,&array); PCHKERRQ(x,ierr);
196 ierr = VecGetLocalSize(x,&n); PCHKERRQ(x,ierr);
198 ierr = VecRestoreArrayRead(x,&array); PCHKERRQ(x,ierr);
204 ierr = VecGetSize(x,&N); PCHKERRQ(x,ierr);
208 PetscParVector::PetscParVector(MPI_Comm comm,
const Vector &_x,
211 ierr = VecCreate(comm,&
x); CCHKERRQ(comm,ierr);
212 ierr = VecSetSizes(
x,_x.
Size(),PETSC_DECIDE); PCHKERRQ(
x,ierr);
213 ierr = VecSetType(
x,VECSTANDARD); PCHKERRQ(
x,ierr);
219 ierr = VecGetArray(
x,&array); PCHKERRQ(
x,ierr);
220 for (
int i = 0; i < _x.
Size(); i++) { array[i] = _x[i]; }
221 ierr = VecRestoreArray(
x,&array); PCHKERRQ(
x,ierr);
228 ierr = VecCreate(comm,&
x); CCHKERRQ(comm,ierr);
232 MPI_Comm_rank(comm, &myid);
233 ierr = VecSetSizes(
x,col[myid+1]-col[myid],PETSC_DECIDE); PCHKERRQ(
x,ierr);
237 ierr = VecSetSizes(
x,PETSC_DECIDE,glob_size); PCHKERRQ(
x,ierr);
239 ierr = VecSetType(
x,VECSTANDARD); PCHKERRQ(
x,ierr);
246 ierr = VecDestroy(&x); CCHKERRQ(comm,ierr);
252 MFEM_VERIFY(col,
"Missing distribution");
254 MPI_Comm_rank(comm, &myid);
255 ierr = VecCreateMPIWithArray(comm,1,col[myid+1]-col[myid],glob_size,_data,
256 &
x); CCHKERRQ(comm,ierr)
262 ierr = VecDuplicate(y.
x,&
x); PCHKERRQ(
x,ierr);
267 bool transpose,
bool allocate) :
Vector()
272 ierr = VecCreate(comm,&
x);
274 ierr = VecSetSizes(
x,loc,PETSC_DECIDE);
276 ierr = VecSetType(
x,VECSTANDARD);
283 ierr = VecCreateMPIWithArray(comm,1,loc,PETSC_DECIDE,NULL,
284 &
x); CCHKERRQ(comm,ierr);
290 bool transpose,
bool allocate) :
Vector()
295 ierr = MatCreateVecs(pA,&
x,NULL); PCHKERRQ(pA,ierr);
299 ierr = MatCreateVecs(pA,NULL,&
x); PCHKERRQ(pA,ierr);
303 ierr = VecReplaceArray(
x,NULL); PCHKERRQ(
x,ierr);
312 ierr = PetscObjectReference((
PetscObject)y); PCHKERRQ(y,ierr);
322 MPI_Comm comm = pfes->
GetComm();
323 ierr = VecCreate(comm,&
x); CCHKERRQ(comm,ierr);
325 PetscMPIInt myid = 0;
326 if (!HYPRE_AssumedPartitionCheck())
328 MPI_Comm_rank(comm,&myid);
330 ierr = VecSetSizes(
x,offsets[myid+1]-offsets[myid],PETSC_DECIDE);
332 ierr = VecSetType(
x,VECSTANDARD); PCHKERRQ(
x,ierr);
348 ierr = VecScatterCreateToAll(
x,&scctx,&vout); PCHKERRQ(
x,ierr);
349 ierr = VecScatterBegin(scctx,
x,vout,INSERT_VALUES,SCATTER_FORWARD);
351 ierr = VecScatterEnd(scctx,
x,vout,INSERT_VALUES,SCATTER_FORWARD);
353 ierr = VecScatterDestroy(&scctx); PCHKERRQ(
x,ierr);
354 ierr = VecGetArray(vout,&array); PCHKERRQ(
x,ierr);
355 ierr = VecGetLocalSize(vout,&size); PCHKERRQ(
x,ierr);
358 ierr = VecRestoreArray(vout,&array); PCHKERRQ(
x,ierr);
359 ierr = VecDestroy(&vout); PCHKERRQ(
x,ierr);
368 ierr = VecSet(
x,d); PCHKERRQ(
x,ierr);
375 MFEM_VERIFY(idx.
Size() == vals.
Size(),
376 "Size mismatch between indices and values");
380 ierr = VecAssemblyBegin(
x); PCHKERRQ(
x,ierr);
381 ierr = VecAssemblyEnd(
x); PCHKERRQ(
x,ierr);
388 MFEM_VERIFY(idx.
Size() == vals.
Size(),
389 "Size mismatch between indices and values");
393 ierr = VecAssemblyBegin(
x); PCHKERRQ(
x,ierr);
394 ierr = VecAssemblyEnd(
x); PCHKERRQ(
x,ierr);
400 ierr = VecCopy(y.
x,
x); PCHKERRQ(
x,ierr);
406 ierr = VecAXPY(
x,1.0,y.
x); PCHKERRQ(
x,ierr);
412 ierr = VecAXPY(
x,-1.0,y.
x); PCHKERRQ(
x,ierr);
418 ierr = VecScale(
x,s); PCHKERRQ(
x,ierr);
424 ierr = VecShift(
x,s); PCHKERRQ(
x,ierr);
430 ierr = VecPlaceArray(
x,temp_data); PCHKERRQ(
x,ierr);
435 ierr = VecResetArray(
x); PCHKERRQ(
x,ierr);
440 PetscRandom rctx = NULL;
444 ierr = PetscRandomCreate(PetscObjectComm((
PetscObject)
x),&rctx);
446 ierr = PetscRandomSetSeed(rctx,(
unsigned long)seed); PCHKERRQ(x,ierr);
447 ierr = PetscRandomSeed(rctx); PCHKERRQ(x,ierr);
449 ierr = VecSetRandom(
x,rctx); PCHKERRQ(
x,ierr);
450 ierr = PetscRandomDestroy(&rctx); PCHKERRQ(
x,ierr);
461 ierr = PetscViewerBinaryOpen(PetscObjectComm((
PetscObject)
x),fname,
462 FILE_MODE_WRITE,&view);
466 ierr = PetscViewerASCIIOpen(PetscObjectComm((
PetscObject)
x),fname,&view);
469 ierr = VecView(
x,view); PCHKERRQ(
x,ierr);
470 ierr = PetscViewerDestroy(&view); PCHKERRQ(
x,ierr);
474 ierr = VecView(
x,NULL); PCHKERRQ(
x,ierr);
483 ierr = MatGetOwnershipRange(
A,&N,NULL); PCHKERRQ(
A,ierr);
490 ierr = MatGetOwnershipRangeColumn(
A,&N,NULL); PCHKERRQ(
A,ierr);
497 ierr = MatGetLocalSize(
A,&N,NULL); PCHKERRQ(
A,ierr);
504 ierr = MatGetLocalSize(
A,NULL,&N); PCHKERRQ(
A,ierr);
511 ierr = MatGetSize(
A,&N,NULL); PCHKERRQ(
A,ierr);
518 ierr = MatGetSize(
A,NULL,&N); PCHKERRQ(
A,ierr);
525 ierr = MatGetInfo(
A,MAT_GLOBAL_SUM,&info); PCHKERRQ(
A,ierr);
550 rows.
GetData(),PETSC_USE_POINTER,&isr); PCHKERRQ(B,ierr);
552 cols.
GetData(),PETSC_USE_POINTER,&isc); PCHKERRQ(B,ierr);
553 ierr = MatCreateSubMatrix(B,isr,isc,MAT_INITIAL_MATRIX,&
A); PCHKERRQ(B,ierr);
554 ierr = ISDestroy(&isr); PCHKERRQ(B,ierr);
555 ierr = ISDestroy(&isc); PCHKERRQ(B,ierr);
599 BlockDiagonalConstructor(comm,row_starts,row_starts,diag,
612 BlockDiagonalConstructor(comm,row_starts,col_starts,diag,
624 ierr = MatDestroy(&A); CCHKERRQ(comm,ierr);
631 #if defined(PETSC_HAVE_HYPRE)
632 ierr = MatCreateFromParCSR(B,MATAIJ,PETSC_USE_POINTER,&
A);
634 ierr = MatConvert_hypreParCSR_AIJ(B,&
A); CCHKERRQ(B.
GetComm(),ierr);
644 ierr = MatDestroy(&A); CCHKERRQ(comm,ierr);
651 ierr = MatDuplicate(B,MAT_COPY_VALUES,&
A); CCHKERRQ(B.
GetComm(),ierr);
659 ierr = MatDuplicate(B,MAT_COPY_VALUES,&
A); CCHKERRQ(B.
GetComm(),ierr);
663 MFEM_VERIFY(
height == B.
Height(),
"Invalid number of local rows");
664 MFEM_VERIFY(
width == B.
Width(),
"Invalid number of local columns");
665 ierr = MatAXPY(
A,1.0,B,DIFFERENT_NONZERO_PATTERN); CCHKERRQ(B.
GetComm(),ierr);
674 ierr = MatDuplicate(B,MAT_COPY_VALUES,&
A); CCHKERRQ(B.
GetComm(),ierr);
675 ierr = MatScale(
A,-1.0); PCHKERRQ(
A,ierr);
679 MFEM_VERIFY(
height == B.
Height(),
"Invalid number of local rows");
680 MFEM_VERIFY(
width == B.
Width(),
"Invalid number of local columns");
681 ierr = MatAXPY(
A,-1.0,B,DIFFERENT_NONZERO_PATTERN); CCHKERRQ(B.
GetComm(),ierr);
686 void PetscParMatrix::
687 BlockDiagonalConstructor(MPI_Comm comm,
692 PetscInt lrsize,lcsize,rstart,cstart;
693 PetscMPIInt myid = 0,commsize;
695 ierr = MPI_Comm_size(comm,&commsize); CCHKERRQ(comm,ierr);
696 if (!HYPRE_AssumedPartitionCheck())
698 ierr = MPI_Comm_rank(comm,&myid); CCHKERRQ(comm,ierr);
700 lrsize = row_starts[myid+1]-row_starts[myid];
701 rstart = row_starts[myid];
702 lcsize = col_starts[myid+1]-col_starts[myid];
703 cstart = col_starts[myid];
708 ierr = ISCreateStride(comm,diag->
Height(),rstart,1,&is); CCHKERRQ(comm,ierr);
709 ISLocalToGlobalMapping rl2g,cl2g;
710 ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g); PCHKERRQ(is,ierr);
711 ierr = ISDestroy(&is); CCHKERRQ(comm,ierr);
712 if (row_starts != col_starts)
714 ierr = ISCreateStride(comm,diag->
Width(),cstart,1,&is);
716 ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g); PCHKERRQ(is,ierr);
717 ierr = ISDestroy(&is); CCHKERRQ(comm,ierr);
721 ierr = PetscObjectReference((
PetscObject)rl2g); PCHKERRQ(rl2g,ierr);
726 ierr = MatCreate(comm,&A); CCHKERRQ(comm,ierr);
727 ierr = MatSetSizes(A,lrsize,lcsize,PETSC_DECIDE,PETSC_DECIDE);
729 ierr = MatSetType(A,MATIS); PCHKERRQ(A,ierr);
730 ierr = MatSetLocalToGlobalMapping(A,rl2g,cl2g); PCHKERRQ(A,ierr);
731 ierr = ISLocalToGlobalMappingDestroy(&rl2g); PCHKERRQ(A,ierr)
732 ierr = ISLocalToGlobalMappingDestroy(&cl2g); PCHKERRQ(A,ierr)
736 ierr = MatISGetLocalMat(A,&lA); PCHKERRQ(A,ierr);
737 int *II = diag->GetI();
738 int *JJ = diag->GetJ();
739 #if defined(PETSC_USE_64BIT_INDICES)
742 ierr = PetscMalloc2(m,&pII,nnz,&pJJ); PCHKERRQ(lA,ierr);
743 for (
int i = 0; i < m; i++) { pII[i] = II[i]; }
744 for (
int i = 0; i < nnz; i++) { pJJ[i] = JJ[i]; }
745 ierr = MatSeqAIJSetPreallocationCSR(lA,pII,pJJ,
746 diag->
GetData()); PCHKERRQ(lA,ierr);
747 ierr = PetscFree2(pII,pJJ); PCHKERRQ(lA,ierr);
749 ierr = MatSeqAIJSetPreallocationCSR(lA,II,JJ,
750 diag->
GetData()); PCHKERRQ(lA,ierr);
762 ierr = PetscMalloc1(m,&dii); CCHKERRQ(PETSC_COMM_SELF,ierr);
763 ierr = PetscMalloc1(nnz,&djj); CCHKERRQ(PETSC_COMM_SELF,ierr);
764 ierr = PetscMalloc1(nnz,&da); CCHKERRQ(PETSC_COMM_SELF,ierr);
765 if (
sizeof(
PetscInt) ==
sizeof(
int))
767 ierr = PetscMemcpy(dii,diag->
GetI(),m*
sizeof(
PetscInt));
768 CCHKERRQ(PETSC_COMM_SELF,ierr);
769 ierr = PetscMemcpy(djj,diag->
GetJ(),nnz*
sizeof(
PetscInt));
770 CCHKERRQ(PETSC_COMM_SELF,ierr);
774 int *iii = diag->
GetI();
775 int *jjj = diag->
GetJ();
776 for (
int i = 0; i < m; i++) { dii[i] = iii[i]; }
777 for (
int i = 0; i < nnz; i++) { djj[i] = jjj[i]; }
780 CCHKERRQ(PETSC_COMM_SELF,ierr);
781 ierr = PetscCalloc1(m,&oii);
782 CCHKERRQ(PETSC_COMM_SELF,ierr);
785 ierr = MatCreateMPIAIJWithSplitArrays(comm,lrsize,lcsize,PETSC_DECIDE,
787 dii,djj,da,oii,NULL,NULL,&A);
792 ierr = MatCreateSeqAIJWithArrays(comm,lrsize,lcsize,dii,djj,da,&A);
796 void *ptrs[4] = {dii,djj,da,oii};
797 const char *names[4] = {
"_mfem_csr_dii",
806 ierr = PetscContainerCreate(comm,&c); CCHKERRQ(comm,ierr);
807 ierr = PetscContainerSetPointer(c,ptrs[i]); CCHKERRQ(comm,ierr);
808 ierr = PetscContainerSetUserDestroy(c,__mfem_array_container_destroy);
812 ierr = PetscContainerDestroy(&c); CCHKERRQ(comm,ierr);
817 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); PCHKERRQ(A,ierr);
818 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); PCHKERRQ(A,ierr);
825 return A ? PetscObjectComm((
PetscObject)A) : MPI_COMM_NULL;
838 ierr = MatCreate(comm,A); CCHKERRQ(comm,ierr);
840 PETSC_DECIDE,PETSC_DECIDE); PCHKERRQ(A,ierr);
841 ierr = MatSetType(*A,MATSHELL); PCHKERRQ(A,ierr);
842 ierr = MatShellSetContext(*A,(
void *)op); PCHKERRQ(A,ierr);
843 ierr = MatShellSetOperation(*A,MATOP_MULT,
844 (
void (*)())__mfem_mat_shell_apply);
846 ierr = MatShellSetOperation(*A,MATOP_MULT_TRANSPOSE,
847 (
void (*)())__mfem_mat_shell_apply_transpose);
849 ierr = MatShellSetOperation(*A,MATOP_COPY,
850 (
void (*)())__mfem_mat_shell_copy);
852 ierr = MatShellSetOperation(*A,MATOP_DESTROY,
853 (
void (*)())__mfem_mat_shell_destroy);
855 ierr = MatSetUp(*A); PCHKERRQ(*A,ierr);
880 PetscBool avoidmatconvert = PETSC_FALSE;
883 ierr = PetscObjectTypeCompareAny((
PetscObject)(pA->
A),&avoidmatconvert,MATMFFD,
887 if (pA && !avoidmatconvert)
891 #if PETSC_VERSION_LT(3,10,0)
895 ierr = PetscObjectTypeCompare((
PetscObject)(pA->
A),MATTRANSPOSEMAT,&istrans);
906 #if PETSC_VERSION_LT(3,10,0)
907 ierr = PetscObjectTypeCompare((
PetscObject)(pA->
A),MATIS,&ismatis);
913 ierr = MatTransposeGetMat(pA->
A,&At); CCHKERRQ(pA->
GetComm(),ierr);
914 #if PETSC_VERSION_LT(3,10,0)
915 ierr = PetscObjectTypeCompare((
PetscObject)(At),MATIS,&ismatis);
923 #if PETSC_VERSION_LT(3,10,0)
930 ierr = MatISGetMPIXAIJ(At,MAT_INITIAL_MATRIX,&B); PCHKERRQ(pA->
A,ierr);
931 ierr = MatCreateTranspose(B,A); PCHKERRQ(pA->
A,ierr);
932 ierr = MatDestroy(&B); PCHKERRQ(pA->
A,ierr);
936 ierr = MatISGetMPIXAIJ(pA->
A,MAT_INITIAL_MATRIX,A);
937 PCHKERRQ(pA->
A,ierr);
944 ierr = MPI_Comm_size(comm,&size); CCHKERRQ(comm,ierr);
950 ierr = MatConvert(At,size > 1 ? MATMPIAIJ : MATSEQAIJ,MAT_INITIAL_MATRIX,&B);
951 PCHKERRQ(pA->
A,ierr);
952 ierr = MatCreateTranspose(B,A); PCHKERRQ(pA->
A,ierr);
953 ierr = MatDestroy(&B); PCHKERRQ(pA->
A,ierr);
957 ierr = MatConvert(pA->
A, size > 1 ? MATMPIAIJ : MATSEQAIJ,MAT_INITIAL_MATRIX,A);
958 PCHKERRQ(pA->
A,ierr);
967 ierr = MatConvert(At,MATIS,MAT_INITIAL_MATRIX,&B); PCHKERRQ(pA->
A,ierr);
968 ierr = MatCreateTranspose(B,A); PCHKERRQ(pA->
A,ierr);
969 ierr = MatDestroy(&B); PCHKERRQ(pA->
A,ierr);
973 ierr = MatConvert(pA->
A,MATIS,MAT_INITIAL_MATRIX,A); PCHKERRQ(pA->
A,ierr);
978 #if defined(PETSC_HAVE_HYPRE)
982 ierr = MatConvert(At,MATHYPRE,MAT_INITIAL_MATRIX,&B); PCHKERRQ(pA->
A,ierr);
983 ierr = MatCreateTranspose(B,A); PCHKERRQ(pA->
A,ierr);
984 ierr = MatDestroy(&B); PCHKERRQ(pA->
A,ierr);
988 ierr = MatConvert(pA->
A,MATHYPRE,MAT_INITIAL_MATRIX,A); PCHKERRQ(pA->
A,ierr);
991 MFEM_ABORT(
"Reconfigure PETSc with --download-hypre or --with-hypre")
1000 MFEM_ABORT(
"Unsupported operator type conversion " << tid)
1007 #if defined(PETSC_HAVE_HYPRE)
1008 ierr = MatCreateFromParCSR(const_cast<HypreParMatrix&>(*pH),MATAIJ,
1009 PETSC_USE_POINTER,A);
1011 ierr = MatConvert_hypreParCSR_AIJ(const_cast<HypreParMatrix&>(*pH),A);
1017 #if defined(PETSC_HAVE_HYPRE)
1018 ierr = MatCreateFromParCSR(const_cast<HypreParMatrix&>(*pH),MATIS,
1019 PETSC_USE_POINTER,A);
1021 ierr = MatConvert_hypreParCSR_IS(const_cast<HypreParMatrix&>(*pH),A);
1027 #if defined(PETSC_HAVE_HYPRE)
1028 ierr = MatCreateFromParCSR(const_cast<HypreParMatrix&>(*pH),MATHYPRE,
1029 PETSC_USE_POINTER,A);
1032 MFEM_ABORT(
"Reconfigure PETSc with --download-hypre or --with-hypre")
1041 MFEM_ABORT(
"Conversion from HypreParCSR to operator type = " << tid <<
1042 " is not implemented");
1047 Mat *mats,*matsl2l = NULL;
1052 ierr = PetscCalloc1(nr*nc,&mats); CCHKERRQ(PETSC_COMM_SELF,ierr);
1055 ierr = PetscCalloc1(nr,&matsl2l); CCHKERRQ(PETSC_COMM_SELF,ierr);
1057 for (i=0; i<nr; i++)
1059 PetscBool needl2l = PETSC_TRUE;
1061 for (j=0; j<nc; j++)
1069 ierr = PetscObjectQuery((
PetscObject)mats[i*nc+j],
"_MatIS_PtAP_l2l",
1071 PCHKERRQ(mats[i*nc+j],ierr);
1079 ierr = PetscContainerGetPointer(c,(
void**)&l2l);
1081 MFEM_VERIFY(l2l->
Size() == 1,
"Unexpected size "
1082 << l2l->
Size() <<
" for block row " << i );
1083 ierr = PetscObjectReference((
PetscObject)(*l2l)[0]);
1085 matsl2l[i] = (*l2l)[0];
1086 needl2l = PETSC_FALSE;
1092 ierr = MatCreateNest(comm,nr,NULL,nc,NULL,mats,A); CCHKERRQ(comm,ierr);
1095 ierr = MatConvert(*A,MATIS,MAT_INPLACE_MATRIX,A); CCHKERRQ(comm,ierr);
1098 for (
int i=0; i<(int)nr; i++) { (*vmatsl2l)[i] = matsl2l[i]; }
1099 ierr = PetscFree(matsl2l); CCHKERRQ(PETSC_COMM_SELF,ierr);
1102 ierr = PetscContainerCreate(comm,&c); CCHKERRQ(comm,ierr);
1103 ierr = PetscContainerSetPointer(c,vmatsl2l); PCHKERRQ(c,ierr);
1104 ierr = PetscContainerSetUserDestroy(c,__mfem_matarray_container_destroy);
1107 PCHKERRQ((*A),ierr);
1108 ierr = PetscContainerDestroy(&c); CCHKERRQ(comm,ierr);
1110 for (i=0; i<nr*nc; i++) { ierr = MatDestroy(&mats[i]); CCHKERRQ(comm,ierr); }
1111 ierr = PetscFree(mats); CCHKERRQ(PETSC_COMM_SELF,ierr);
1117 ierr = MatCreate(comm,A); CCHKERRQ(comm,ierr);
1118 ierr = MatSetSizes(*A,pI->
Height(),pI->
Width(),PETSC_DECIDE,PETSC_DECIDE);
1120 ierr = MatSetType(*A,MATAIJ); PCHKERRQ(*A,ierr);
1121 ierr = MatMPIAIJSetPreallocation(*A,1,NULL,0,NULL); PCHKERRQ(*A,ierr);
1122 ierr = MatSeqAIJSetPreallocation(*A,1,NULL); PCHKERRQ(*A,ierr);
1123 ierr = MatSetOption(*A,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE); PCHKERRQ(*A,ierr);
1124 ierr = MatGetOwnershipRange(*A,&rst,NULL); PCHKERRQ(*A,ierr);
1127 ierr = MatSetValue(*A,i,i,1.,INSERT_VALUES); PCHKERRQ(*A,ierr);
1129 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY); PCHKERRQ(*A,ierr);
1130 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY); PCHKERRQ(*A,ierr);
1146 int n = pS->
Width();
1147 int *ii = pS->
GetI();
1148 int *jj = pS->
GetJ();
1151 ierr = PetscMalloc1(m+1,&pii); CCHKERRQ(PETSC_COMM_SELF,ierr);
1152 ierr = PetscMalloc1(ii[m],&pjj); CCHKERRQ(PETSC_COMM_SELF,ierr);
1153 ierr = PetscMalloc1(ii[m],&pdata); CCHKERRQ(PETSC_COMM_SELF,ierr);
1155 for (
int i = 0; i < m; i++)
1157 bool issorted =
true;
1159 for (
int j = ii[i]; j < ii[i+1]; j++)
1162 if (j != ii[i] && issorted) { issorted = (pjj[j] > pjj[j-1]); }
1167 ierr = PetscSortIntWithScalarArray(pii[i+1]-pii[i],pjj + pii[i],pdata + pii[i]);
1168 CCHKERRQ(PETSC_COMM_SELF,ierr);
1172 ierr = MPI_Comm_size(comm,&size); CCHKERRQ(comm,ierr);
1175 ierr = MatCreateSeqAIJWithArrays(comm,m,n,pii,pjj,pdata,&B);
1176 CCHKERRQ(comm,ierr);
1181 ierr = PetscCalloc1(m+1,&oii); CCHKERRQ(PETSC_COMM_SELF,ierr);
1182 ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,
1184 pii,pjj,pdata,oii,NULL,NULL,&B);
1185 CCHKERRQ(comm,ierr);
1187 void *ptrs[4] = {pii,pjj,pdata,oii};
1188 const char *names[4] = {
"_mfem_csr_pii",
1193 for (
int i=0; i<4; i++)
1197 ierr = PetscContainerCreate(PETSC_COMM_SELF,&c); PCHKERRQ(B,ierr);
1198 ierr = PetscContainerSetPointer(c,ptrs[i]); PCHKERRQ(B,ierr);
1199 ierr = PetscContainerSetUserDestroy(c,__mfem_array_container_destroy);
1203 ierr = PetscContainerDestroy(&c); PCHKERRQ(B,ierr);
1211 ierr = MatConvert(B,MATHYPRE,MAT_INITIAL_MATRIX,A); PCHKERRQ(B,ierr);
1212 ierr = MatDestroy(&B); PCHKERRQ(*A,ierr);
1216 ierr = MatConvert(B,MATIS,MAT_INITIAL_MATRIX,A); PCHKERRQ(B,ierr);
1217 ierr = MatDestroy(&B); PCHKERRQ(*A,ierr);
1221 MFEM_ABORT(
"Unsupported operator type conversion " << tid)
1228 "Supported types are ANY_TYPE, PETSC_MATSHELL or PETSC_MATAIJ");
1235 ierr = MatComputeOperator(*A,MATMPIAIJ,&B); CCHKERRQ(comm,ierr);
1236 ierr = PetscObjectTypeCompare((
PetscObject)B,MATMPIAIJ,&isaij);
1237 CCHKERRQ(comm,ierr);
1238 ierr = MatDestroy(A); CCHKERRQ(comm,ierr);
1241 ierr = MatConvert(B,MATAIJ,MAT_INITIAL_MATRIX,A); CCHKERRQ(comm,ierr);
1242 ierr = MatDestroy(&B); CCHKERRQ(comm,ierr);
1256 MPI_Comm comm = MPI_COMM_NULL;
1257 ierr = PetscObjectGetComm((
PetscObject)A,&comm); PCHKERRQ(A,ierr);
1258 ierr = MatDestroy(&A); CCHKERRQ(comm,ierr);
1269 ierr = PetscObjectReference((
PetscObject)a); PCHKERRQ(a,ierr);
1279 if (_A == A) {
return; }
1281 ierr = PetscObjectReference((
PetscObject)_A); PCHKERRQ(_A,ierr);
1296 f = MatMultTranspose;
1297 fadd = MatMultTransposeAdd;
1308 ierr = VecScale(X,a); PCHKERRQ(A,ierr);
1309 ierr = (*fadd)(A,X,Y,Y); PCHKERRQ(A,ierr);
1310 ierr = VecScale(X,1./a); PCHKERRQ(A,ierr);
1314 ierr = VecScale(X,a); PCHKERRQ(A,ierr);
1315 ierr = VecScale(Y,b); PCHKERRQ(A,ierr);
1316 ierr = (*fadd)(A,X,Y,Y); PCHKERRQ(A,ierr);
1317 ierr = VecScale(X,1./a); PCHKERRQ(A,ierr);
1321 ierr = (*f)(A,X,Y); PCHKERRQ(A,ierr);
1324 ierr = VecScale(Y,a); PCHKERRQ(A,ierr);
1336 ierr = VecScale(Y,b); PCHKERRQ(A,ierr);
1340 ierr = VecSet(Y,0.); PCHKERRQ(A,ierr);
1347 ierr = PetscObjectReference((
PetscObject)master.
A); PCHKERRQ(master.
A,ierr);
1359 MFEM_VERIFY(A,
"Mat not present");
1369 MFEM_VERIFY(A,
"Mat not present");
1380 ierr = MatCreateTranspose(A,&B); PCHKERRQ(A,ierr);
1384 ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&B); PCHKERRQ(A,ierr);
1391 ierr = MatScale(A,s); PCHKERRQ(A,ierr);
1396 MFEM_ASSERT(x.
Size() ==
Width(),
"invalid x.Size() = " << x.
Size()
1397 <<
", expected size = " <<
Width());
1398 MFEM_ASSERT(y.
Size() ==
Height(),
"invalid y.Size() = " << y.
Size()
1399 <<
", expected size = " <<
Height());
1405 MatMultKernel(A,a,XX->x,b,YY->
x,
false);
1413 MFEM_ASSERT(x.
Size() ==
Height(),
"invalid x.Size() = " << x.
Size()
1414 <<
", expected size = " <<
Height());
1415 MFEM_ASSERT(y.
Size() ==
Width(),
"invalid y.Size() = " << y.
Size()
1416 <<
", expected size = " <<
Width());
1422 MatMultKernel(A,a,YY->
x,b,XX->x,
true);
1435 ierr = PetscViewerBinaryOpen(PetscObjectComm((
PetscObject)A),fname,
1436 FILE_MODE_WRITE,&view);
1440 ierr = PetscViewerASCIIOpen(PetscObjectComm((
PetscObject)A),fname,&view);
1443 ierr = MatView(A,view); PCHKERRQ(A,ierr);
1444 ierr = PetscViewerDestroy(&view); PCHKERRQ(A,ierr);
1448 ierr = MatView(A,NULL); PCHKERRQ(A,ierr);
1454 MFEM_ASSERT(s.
Size() ==
Height(),
"invalid s.Size() = " << s.
Size()
1455 <<
", expected size = " <<
Height());
1459 ierr = MatDiagonalScale(A,*YY,NULL); PCHKERRQ(A,ierr);
1465 MFEM_ASSERT(s.
Size() ==
Width(),
"invalid s.Size() = " << s.
Size()
1466 <<
", expected size = " <<
Width());
1470 ierr = MatDiagonalScale(A,NULL,*XX); PCHKERRQ(A,ierr);
1476 ierr = MatShift(A,(
PetscScalar)s); PCHKERRQ(A,ierr);
1482 MFEM_ASSERT(s.
Size() ==
Height(),
"invalid s.Size() = " << s.
Size()
1483 <<
", expected size = " <<
Height());
1484 MFEM_ASSERT(s.
Size() ==
Width(),
"invalid s.Size() = " << s.
Size()
1485 <<
", expected size = " <<
Width());
1489 ierr = MatDiagonalSet(A,*XX,ADD_VALUES); PCHKERRQ(A,ierr);
1497 "Petsc TripleMatrixProduct: Number of local cols of A " << A->
Width() <<
1498 " differs from number of local rows of P " << P->
Height());
1500 "Petsc TripleMatrixProduct: Number of local rows of A " << A->
Height() <<
1501 " differs from number of local cols of R " << R->
Width());
1503 ierr = MatMatMatMult(*R,*A,*P,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B);
1510 Mat pA = *A,pP = *P,pRt = *Rt;
1512 PetscBool Aismatis,Pismatis,Rtismatis;
1515 "Petsc RAP: Number of local cols of A " << A->
Width() <<
1516 " differs from number of local rows of P " << P->
Height());
1518 "Petsc RAP: Number of local rows of A " << A->
Height() <<
1519 " differs from number of local rows of Rt " << Rt->
Height());
1520 ierr = PetscObjectTypeCompare((
PetscObject)pA,MATIS,&Aismatis);
1522 ierr = PetscObjectTypeCompare((
PetscObject)pP,MATIS,&Pismatis);
1524 ierr = PetscObjectTypeCompare((
PetscObject)pRt,MATIS,&Rtismatis);
1531 ISLocalToGlobalMapping cl2gP,cl2gRt;
1532 PetscInt rlsize,clsize,rsize,csize;
1534 ierr = MatGetLocalToGlobalMapping(pP,NULL,&cl2gP); PCHKERRQ(pA,ierr);
1535 ierr = MatGetLocalToGlobalMapping(pRt,NULL,&cl2gRt); PCHKERRQ(pA,ierr);
1536 ierr = MatGetLocalSize(pP,NULL,&clsize); PCHKERRQ(pP,ierr);
1537 ierr = MatGetLocalSize(pRt,NULL,&rlsize); PCHKERRQ(pRt,ierr);
1538 ierr = MatGetSize(pP,NULL,&csize); PCHKERRQ(pP,ierr);
1539 ierr = MatGetSize(pRt,NULL,&rsize); PCHKERRQ(pRt,ierr);
1540 ierr = MatCreate(A->
GetComm(),&B); PCHKERRQ(pA,ierr);
1541 ierr = MatSetSizes(B,rlsize,clsize,rsize,csize); PCHKERRQ(B,ierr);
1542 ierr = MatSetType(B,MATIS); PCHKERRQ(B,ierr);
1543 ierr = MatSetLocalToGlobalMapping(B,cl2gRt,cl2gP); PCHKERRQ(B,ierr);
1544 ierr = MatISGetLocalMat(pA,&lA); PCHKERRQ(pA,ierr);
1545 ierr = MatISGetLocalMat(pP,&lP); PCHKERRQ(pA,ierr);
1546 ierr = MatISGetLocalMat(pRt,&lRt); PCHKERRQ(pA,ierr);
1549 ierr = MatPtAP(lA,lP,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&lB);
1555 ierr = MatTranspose(lRt,MAT_INITIAL_MATRIX,&lR); PCHKERRQ(lRt,ierr);
1556 ierr = MatMatMatMult(lR,lA,lP,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&lB);
1558 ierr = MatDestroy(&lR); PCHKERRQ(lRt,ierr);
1566 ierr = PetscObjectReference((
PetscObject)lRt); PCHKERRQ(lRt,ierr);
1567 (*vmatsl2l)[0] = lRt;
1570 ierr = PetscContainerCreate(PetscObjectComm((
PetscObject)B),&c);
1572 ierr = PetscContainerSetPointer(c,vmatsl2l); PCHKERRQ(c,ierr);
1573 ierr = PetscContainerSetUserDestroy(c,__mfem_matarray_container_destroy);
1577 ierr = PetscContainerDestroy(&c); PCHKERRQ(B,ierr);
1581 ierr = MatISSetLocalMat(B,lB); PCHKERRQ(lB,ierr);
1582 ierr = MatDestroy(&lB); PCHKERRQ(lA,ierr);
1583 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); PCHKERRQ(B,ierr);
1584 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); PCHKERRQ(B,ierr);
1590 ierr = MatPtAP(pA,pP,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B);
1596 ierr = MatTranspose(pRt,MAT_INITIAL_MATRIX,&pR); PCHKERRQ(Rt,ierr);
1597 ierr = MatMatMatMult(pR,pA,pP,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B);
1599 ierr = MatDestroy(&pR); PCHKERRQ(pRt,ierr);
1614 #if defined(PETSC_HAVE_HYPRE)
1629 ierr = MatMatMult(*A,*B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AB);
1639 ierr = MatDuplicate(A,MAT_COPY_VALUES,&Ae); PCHKERRQ(A,ierr);
1641 ierr = MatAXPY(Ae,-1.,A,SAME_NONZERO_PATTERN); PCHKERRQ(A,ierr);
1650 MFEM_ABORT(
"Missing PetscParMatrix::EliminateRowsCols() with HypreParVectors");
1659 ierr = MatGetSize(A,&M,&N); PCHKERRQ(A,ierr);
1660 MFEM_VERIFY(M == N,
"Rectangular case unsupported");
1663 ierr = MatSetOption(A,MAT_NO_OFF_PROC_ZERO_ROWS,PETSC_TRUE); PCHKERRQ(A,ierr);
1667 ierr = MatGetOwnershipRange(A,&rst,NULL); PCHKERRQ(A,ierr);
1670 ierr = Convert_Array_IS(
GetComm(),
true,&rows_cols,rst,&dir); PCHKERRQ(A,ierr);
1673 ierr = MatZeroRowsColumnsIS(A,dir,diag,NULL,NULL); PCHKERRQ(A,ierr);
1677 ierr = MatZeroRowsColumnsIS(A,dir,diag,X,B); PCHKERRQ(A,ierr);
1679 ierr = ISDestroy(&dir); PCHKERRQ(A,ierr);
1684 ierr = MatSetOption(A,MAT_NO_OFF_PROC_ZERO_ROWS,PETSC_TRUE); PCHKERRQ(A,ierr);
1688 ierr = MatGetOwnershipRange(A,&rst,NULL); PCHKERRQ(A,ierr);
1691 ierr = Convert_Array_IS(
GetComm(),
true,&rows,rst,&dir); PCHKERRQ(A,ierr);
1692 ierr = MatZeroRowsIS(A,dir,0.0,NULL,NULL); PCHKERRQ(A,ierr);
1693 ierr = ISDestroy(&dir); PCHKERRQ(A,ierr);
1703 ierr = PetscObjectDereference((
PetscObject)A); CCHKERRQ(comm,ierr);
1712 MFEM_VERIFY(A,
"no associated PETSc Mat object");
1715 ierr = PetscObjectTypeCompare(oA, MATAIJ, &ok); PCHKERRQ(A,ierr);
1717 ierr = PetscObjectTypeCompare(oA, MATSEQAIJ, &ok); PCHKERRQ(A,ierr);
1719 ierr = PetscObjectTypeCompare(oA, MATMPIAIJ, &ok); PCHKERRQ(A,ierr);
1721 ierr = PetscObjectTypeCompare(oA, MATIS, &ok); PCHKERRQ(A,ierr);
1723 ierr = PetscObjectTypeCompare(oA, MATSHELL, &ok); PCHKERRQ(A,ierr);
1725 ierr = PetscObjectTypeCompare(oA, MATNEST, &ok); PCHKERRQ(A,ierr);
1727 #if defined(PETSC_HAVE_HYPRE)
1728 ierr = PetscObjectTypeCompare(oA, MATHYPRE, &ok); PCHKERRQ(A,ierr);
1742 Ae.
Mult(-1.0, X, 1.0, B);
1745 ierr = MatGetDiagonal(pA,diag); PCHKERRQ(pA,ierr);
1746 ierr = VecGetArrayRead(diag,&array); PCHKERRQ(diag,ierr);
1747 for (
int i = 0; i < ess_dof_list.
Size(); i++)
1749 int r = ess_dof_list[i];
1750 B(r) = array[r] * X(r);
1752 ierr = VecRestoreArrayRead(diag,&array); PCHKERRQ(diag,ierr);
1781 if (
cid == KSP_CLASSID)
1784 ierr = KSPSetTolerances(ksp,tol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
1786 else if (
cid == SNES_CLASSID)
1789 ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tol,PETSC_DEFAULT,PETSC_DEFAULT,
1792 else if (
cid == TS_CLASSID)
1795 ierr = TSSetTolerances(ts,PETSC_DECIDE,NULL,tol,NULL);
1799 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
1806 if (
cid == KSP_CLASSID)
1809 ierr = KSPSetTolerances(ksp,PETSC_DEFAULT,tol,PETSC_DEFAULT,PETSC_DEFAULT);
1811 else if (
cid == SNES_CLASSID)
1814 ierr = SNESSetTolerances(snes,tol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,
1817 else if (
cid == TS_CLASSID)
1820 ierr = TSSetTolerances(ts,tol,NULL,PETSC_DECIDE,NULL);
1824 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
1831 if (
cid == KSP_CLASSID)
1834 ierr = KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,
1837 else if (
cid == SNES_CLASSID)
1840 ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,
1841 max_iter,PETSC_DEFAULT);
1843 else if (
cid == TS_CLASSID)
1846 ierr = TSSetMaxSteps(ts,max_iter);
1850 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
1858 typedef PetscErrorCode (*myPetscFunc)(
void**);
1859 PetscViewerAndFormat *vf = NULL;
1860 PetscViewer viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm(
obj));
1864 ierr = PetscViewerAndFormatCreate(viewer,PETSC_VIEWER_DEFAULT,&vf);
1867 if (
cid == KSP_CLASSID)
1875 ierr = KSPMonitorCancel(ksp); PCHKERRQ(ksp,ierr);
1879 ierr = KSPMonitorSet(ksp,(myMonitor)KSPMonitorDefault,vf,
1880 (myPetscFunc)PetscViewerAndFormatDestroy);
1885 ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE); PCHKERRQ(ksp,ierr);
1886 ierr = KSPMonitorSet(ksp,(myMonitor)KSPMonitorSingularValue,vf,
1887 (myPetscFunc)PetscViewerAndFormatDestroy);
1891 ierr = PetscViewerAndFormatCreate(viewer,PETSC_VIEWER_DEFAULT,&vf);
1892 PCHKERRQ(viewer,ierr);
1893 ierr = KSPMonitorSet(ksp,(myMonitor)KSPMonitorTrueResidualNorm,vf,
1894 (myPetscFunc)PetscViewerAndFormatDestroy);
1899 else if (
cid == SNES_CLASSID)
1905 ierr = SNESMonitorCancel(snes); PCHKERRQ(snes,ierr);
1909 ierr = SNESMonitorSet(snes,(myMonitor)SNESMonitorDefault,vf,
1910 (myPetscFunc)PetscViewerAndFormatDestroy);
1911 PCHKERRQ(snes,ierr);
1914 else if (
cid == TS_CLASSID)
1919 ierr = TSMonitorCancel(ts); PCHKERRQ(ts,ierr);
1924 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
1930 return obj ? PetscObjectComm(
obj) : MPI_COMM_NULL;
1935 __mfem_monitor_ctx *monctx;
1936 ierr = PetscNew(&monctx); CCHKERRQ(PETSC_COMM_SELF,ierr);
1937 monctx->solver =
this;
1938 monctx->monitor = ctx;
1939 if (
cid == KSP_CLASSID)
1941 ierr = KSPMonitorSet((
KSP)
obj,__mfem_ksp_monitor,monctx,
1942 __mfem_monitor_ctx_destroy);
1945 else if (
cid == SNES_CLASSID)
1947 ierr = SNESMonitorSet((
SNES)
obj,__mfem_snes_monitor,monctx,
1948 __mfem_monitor_ctx_destroy);
1951 else if (
cid == TS_CLASSID)
1953 ierr = TSMonitorSet((
TS)
obj,__mfem_ts_monitor,monctx,
1954 __mfem_monitor_ctx_destroy);
1959 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
1966 if (
cid == SNES_CLASSID)
1968 __mfem_snes_ctx* snes_ctx = (__mfem_snes_ctx*)
private_ctx;
1971 else if (
cid == TS_CLASSID)
1973 __mfem_ts_ctx* ts_ctx = (__mfem_ts_ctx*)
private_ctx;
1978 MFEM_ABORT(
"Handling of essential bc only implemented for nonlinear and time-dependent solvers");
1985 if (
cid == TS_CLASSID)
1990 ierr = TSGetSNES((
TS)
obj,&snes); PCHKERRQ(obj,ierr);
1991 ierr = SNESGetKSP(snes,&ksp); PCHKERRQ(obj,ierr);
1992 ierr = KSPGetPC(ksp,&pc); PCHKERRQ(obj,ierr);
1994 else if (
cid == SNES_CLASSID)
1998 ierr = SNESGetKSP((
SNES)
obj,&ksp); PCHKERRQ(obj,ierr);
1999 ierr = KSPGetPC(ksp,&pc); PCHKERRQ(obj,ierr);
2001 else if (
cid == KSP_CLASSID)
2003 ierr = KSPGetPC((
KSP)
obj,&pc); PCHKERRQ(obj,ierr);
2005 else if (
cid == PC_CLASSID)
2011 MFEM_ABORT(
"No support for PetscPreconditionerFactory for this object");
2015 ierr = MakeShellPCWithFactory(pc,factory); PCHKERRQ(pc,ierr);
2019 ierr = PCSetType(pc, PCNONE); PCHKERRQ(pc,ierr);
2025 if (!customize) {
clcustom =
true; }
2028 if (
cid == PC_CLASSID)
2031 ierr = PCSetFromOptions(pc); PCHKERRQ(pc, ierr);
2033 else if (
cid == KSP_CLASSID)
2036 ierr = KSPSetFromOptions(ksp); PCHKERRQ(ksp, ierr);
2038 else if (
cid == SNES_CLASSID)
2041 ierr = SNESSetFromOptions(snes); PCHKERRQ(snes, ierr);
2043 else if (
cid == TS_CLASSID)
2046 ierr = TSSetFromOptions(ts); PCHKERRQ(ts, ierr);
2050 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
2058 if (
cid == KSP_CLASSID)
2061 KSPConvergedReason reason;
2062 ierr = KSPGetConvergedReason(ksp,&reason);
2064 return reason > 0 ? 1 : 0;
2066 else if (
cid == SNES_CLASSID)
2069 SNESConvergedReason reason;
2070 ierr = SNESGetConvergedReason(snes,&reason);
2071 PCHKERRQ(snes,ierr);
2072 return reason > 0 ? 1 : 0;
2074 else if (
cid == TS_CLASSID)
2077 TSConvergedReason reason;
2078 ierr = TSGetConvergedReason(ts,&reason);
2080 return reason > 0 ? 1 : 0;
2084 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
2091 if (
cid == KSP_CLASSID)
2095 ierr = KSPGetIterationNumber(ksp,&its);
2099 else if (
cid == SNES_CLASSID)
2103 ierr = SNESGetIterationNumber(snes,&its);
2104 PCHKERRQ(snes,ierr);
2107 else if (
cid == TS_CLASSID)
2111 ierr = TSGetStepNumber(ts,&its);
2117 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
2124 if (
cid == KSP_CLASSID)
2128 ierr = KSPGetResidualNorm(ksp,&norm);
2132 if (
cid == SNES_CLASSID)
2136 ierr = SNESGetFunctionNorm(snes,&norm);
2137 PCHKERRQ(snes,ierr);
2142 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
2143 return PETSC_MAX_REAL;
2150 if (
cid == SNES_CLASSID)
2152 __mfem_snes_ctx *snes_ctx;
2153 ierr = PetscNew(&snes_ctx); CCHKERRQ(PETSC_COMM_SELF,ierr);
2154 snes_ctx->op = NULL;
2155 snes_ctx->bchandler = NULL;
2156 snes_ctx->work = NULL;
2160 else if (
cid == TS_CLASSID)
2162 __mfem_ts_ctx *ts_ctx;
2163 ierr = PetscNew(&ts_ctx); CCHKERRQ(PETSC_COMM_SELF,ierr);
2165 ts_ctx->bchandler = NULL;
2166 ts_ctx->work = NULL;
2167 ts_ctx->work2 = NULL;
2168 ts_ctx->cached_shift = std::numeric_limits<PetscReal>::min();
2169 ts_ctx->cached_ijacstate = -1;
2170 ts_ctx->cached_rhsjacstate = -1;
2171 ts_ctx->cached_splits_xstate = -1;
2172 ts_ctx->cached_splits_xdotstate = -1;
2183 if (
cid == SNES_CLASSID)
2185 __mfem_snes_ctx *snes_ctx = (__mfem_snes_ctx *)
private_ctx;
2186 delete snes_ctx->work;
2188 else if (
cid == TS_CLASSID)
2190 __mfem_ts_ctx *ts_ctx = (__mfem_ts_ctx *)
private_ctx;
2191 delete ts_ctx->work;
2192 delete ts_ctx->work2;
2194 ierr = PetscFree(
private_ctx); CCHKERRQ(PETSC_COMM_SELF,ierr);
2201 : bctype(_type), setup(false), eval_t(0.0),
2202 eval_t_cached(std::numeric_limits<double>::min())
2210 ess_tdof_list.
Assign(list);
2216 if (setup) {
return; }
2220 this->
Eval(eval_t,eval_g);
2221 eval_t_cached = eval_t;
2232 (*this).SetUp(x.
Size());
2236 for (
int i = 0; i < ess_tdof_list.
Size(); ++i)
2238 y[ess_tdof_list[i]] = 0.0;
2243 if (bctype !=
CONSTANT && eval_t != eval_t_cached)
2245 Eval(eval_t,eval_g);
2246 eval_t_cached = eval_t;
2248 for (
int i = 0; i < ess_tdof_list.
Size(); ++i)
2250 y[ess_tdof_list[i]] = eval_g[ess_tdof_list[i]];
2257 (*this).SetUp(x.
Size());
2260 for (
int i = 0; i < ess_tdof_list.
Size(); ++i)
2262 x[ess_tdof_list[i]] = 0.0;
2267 if (bctype !=
CONSTANT && eval_t != eval_t_cached)
2269 Eval(eval_t,eval_g);
2270 eval_t_cached = eval_t;
2272 for (
int i = 0; i < ess_tdof_list.
Size(); ++i)
2274 x[ess_tdof_list[i]] = eval_g[ess_tdof_list[i]];
2281 (*this).SetUp(x.
Size());
2284 for (
int i = 0; i < ess_tdof_list.
Size(); ++i)
2286 y[ess_tdof_list[i]] = x[ess_tdof_list[i]];
2291 for (
int i = 0; i < ess_tdof_list.
Size(); ++i)
2293 y[ess_tdof_list[i]] = x[ess_tdof_list[i]] - eval_g[ess_tdof_list[i]];
2300 (*this).SetUp(x.
Size());
2301 for (
int i = 0; i < ess_tdof_list.
Size(); ++i)
2303 x[ess_tdof_list[i]] = 0.0;
2309 (*this).SetUp(x.
Size());
2311 for (
int i = 0; i < ess_tdof_list.
Size(); ++i)
2313 y[ess_tdof_list[i]] = 0.0;
2324 ierr = KSPCreate(comm,&ksp); CCHKERRQ(comm,ierr);
2326 ierr = PetscObjectGetClassId(
obj,&
cid); PCHKERRQ(
obj,ierr);
2327 ierr = KSPSetOptionsPrefix(ksp, prefix.c_str()); PCHKERRQ(ksp, ierr);
2331 const std::string &prefix)
2337 ierr = PetscObjectGetClassId(obj,&
cid); PCHKERRQ(obj,ierr);
2338 ierr = KSPSetOptionsPrefix(ksp, prefix.c_str()); PCHKERRQ(ksp, ierr);
2343 const std::string &prefix)
2349 ierr = PetscObjectGetClassId(obj, &
cid); PCHKERRQ(obj, ierr);
2350 ierr = KSPSetOptionsPrefix(ksp, prefix.c_str()); PCHKERRQ(ksp, ierr);
2362 bool delete_pA =
false;
2380 MFEM_VERIFY(pA,
"Unsupported operation!");
2388 PetscInt nheight,nwidth,oheight,owidth;
2390 ierr = KSPGetOperators(ksp,&C,NULL); PCHKERRQ(ksp,ierr);
2391 ierr = MatGetSize(A,&nheight,&nwidth); PCHKERRQ(A,ierr);
2392 ierr = MatGetSize(C,&oheight,&owidth); PCHKERRQ(A,ierr);
2393 if (nheight != oheight || nwidth != owidth)
2397 ierr = KSPReset(ksp); PCHKERRQ(ksp,ierr);
2404 ierr = KSPSetOperators(ksp,A,A); PCHKERRQ(ksp,ierr);
2413 if (delete_pA) {
delete pA; }
2428 bool delete_pA =
false;
2446 MFEM_VERIFY(pA,
"Unsupported operation!");
2449 bool delete_ppA =
false;
2452 if (oA == poA && !wrap)
2462 MFEM_VERIFY(ppA,
"Unsupported operation!");
2471 PetscInt nheight,nwidth,oheight,owidth;
2473 ierr = KSPGetOperators(ksp,&C,NULL); PCHKERRQ(ksp,ierr);
2474 ierr = MatGetSize(A,&nheight,&nwidth); PCHKERRQ(A,ierr);
2475 ierr = MatGetSize(C,&oheight,&owidth); PCHKERRQ(A,ierr);
2476 if (nheight != oheight || nwidth != owidth)
2480 ierr = KSPReset(ksp); PCHKERRQ(ksp,ierr);
2487 ierr = KSPSetOperators(ksp,A,P); PCHKERRQ(ksp,ierr);
2496 if (delete_pA) {
delete pA; }
2497 if (delete_ppA) {
delete ppA; }
2507 ierr = KSPGetOperatorsSet(ksp,&amat,NULL); PCHKERRQ(ksp,ierr);
2510 ierr = KSPGetOperators(ksp,&A,NULL); PCHKERRQ(ksp,ierr);
2511 ierr = PetscObjectReference((
PetscObject)A); PCHKERRQ(ksp,ierr);
2516 ierr = KSPSetPC(ksp,*ppc); PCHKERRQ(ksp,ierr);
2525 ierr = KSPGetPC(ksp,&pc); PCHKERRQ(ksp,ierr);
2526 ierr = MakeShellPC(pc,precond,
false); PCHKERRQ(ksp,ierr);
2532 ierr = KSPGetOperators(ksp,NULL,&P); PCHKERRQ(ksp,ierr);
2533 ierr = PetscObjectReference((
PetscObject)P); PCHKERRQ(ksp,ierr);
2534 ierr = KSPSetOperators(ksp,A,P); PCHKERRQ(ksp,ierr);
2535 ierr = MatDestroy(&A); PCHKERRQ(ksp,ierr);
2536 ierr = MatDestroy(&P); PCHKERRQ(ksp,ierr);
2540 void PetscLinearSolver::MultKernel(
const Vector &b,
Vector &x,
bool trans)
const
2547 ierr = KSPGetOperators(ksp, &pA, NULL); PCHKERRQ(
obj, ierr);
2555 PetscParMatrix A = PetscParMatrix(pA,
true);
2556 X =
new PetscParVector(A,
false,
false);
2564 ierr = KSPSetInitialGuessNonzero(ksp, (PetscBool)
iterative_mode);
2565 PCHKERRQ(ksp, ierr);
2570 ierr = KSPSolveTranspose(ksp,
B->
x, X->x); PCHKERRQ(ksp,ierr);
2574 ierr = KSPSolve(ksp,
B->
x, X->x); PCHKERRQ(ksp,ierr);
2582 (*this).MultKernel(b,x,
false);
2587 (*this).MultKernel(b,x,
true);
2594 ierr = PetscObjectGetComm((
PetscObject)ksp,&comm); PCHKERRQ(ksp,ierr);
2595 ierr = KSPDestroy(&ksp); CCHKERRQ(comm,ierr);
2604 ierr = KSPSetType(ksp,KSPCG); PCHKERRQ(ksp,ierr);
2606 ierr = KSPSetNormType(ksp,KSP_NORM_NATURAL); PCHKERRQ(ksp,ierr);
2613 ierr = KSPSetType(ksp,KSPCG); PCHKERRQ(ksp,ierr);
2615 ierr = KSPSetNormType(ksp,KSP_NORM_NATURAL); PCHKERRQ(ksp,ierr);
2619 const std::string &prefix)
2623 ierr = KSPSetType(ksp,KSPCG); PCHKERRQ(ksp,ierr);
2625 ierr = KSPSetNormType(ksp,KSP_NORM_NATURAL); PCHKERRQ(ksp,ierr);
2631 const std::string &prefix)
2635 ierr = PCCreate(comm,&pc); CCHKERRQ(comm,ierr);
2637 ierr = PetscObjectGetClassId(
obj,&
cid); PCHKERRQ(
obj,ierr);
2638 ierr = PCSetOptionsPrefix(pc, prefix.c_str()); PCHKERRQ(pc, ierr);
2642 const string &prefix)
2648 ierr = PetscObjectGetClassId(obj,&
cid); PCHKERRQ(obj,ierr);
2649 ierr = PCSetOptionsPrefix(pc, prefix.c_str()); PCHKERRQ(pc, ierr);
2654 const string &prefix)
2658 ierr = PCCreate(comm,&pc); CCHKERRQ(comm,ierr);
2660 ierr = PetscObjectGetClassId(
obj,&
cid); PCHKERRQ(
obj,ierr);
2661 ierr = PCSetOptionsPrefix(pc, prefix.c_str()); PCHKERRQ(pc, ierr);
2667 bool delete_pA =
false;
2684 PetscInt nheight,nwidth,oheight,owidth;
2686 ierr = PCGetOperators(pc,&C,NULL); PCHKERRQ(pc,ierr);
2687 ierr = MatGetSize(A,&nheight,&nwidth); PCHKERRQ(A,ierr);
2688 ierr = MatGetSize(C,&oheight,&owidth); PCHKERRQ(A,ierr);
2689 if (nheight != oheight || nwidth != owidth)
2693 ierr = PCReset(pc); PCHKERRQ(pc,ierr);
2699 ierr = PCSetOperators(pc,pA->
A,pA->
A); PCHKERRQ(
obj,ierr);
2708 if (delete_pA) {
delete pA; };
2711 void PetscPreconditioner::MultKernel(
const Vector &b,
Vector &x,
2719 ierr = PCGetOperators(pc, NULL, &pA); PCHKERRQ(
obj, ierr);
2727 PetscParMatrix A(pA,
true);
2728 X =
new PetscParVector(A,
false,
false);
2739 ierr = PCApplyTranspose(pc,
B->
x, X->x); PCHKERRQ(pc, ierr);
2743 ierr = PCApply(pc,
B->
x, X->x); PCHKERRQ(pc, ierr);
2751 (*this).MultKernel(b,x,
false);
2756 (*this).MultKernel(b,x,
true);
2763 ierr = PetscObjectGetComm((
PetscObject)pc,&comm); PCHKERRQ(pc,ierr);
2764 ierr = PCDestroy(&pc); CCHKERRQ(comm,ierr);
2772 for (
int i = 0; i < x.
Size(); i++) { y(i) = x(i); }
2775 void PetscBDDCSolver::BDDCSolverConstructor(
const PetscBDDCSolverParams &opts)
2777 MPI_Comm comm = PetscObjectComm(
obj);
2782 ierr = PCGetOperators(pc,NULL,&pA); PCHKERRQ(pc,ierr);
2786 ierr = PetscObjectTypeCompare((
PetscObject)pA,MATIS,&ismatis);
2788 MFEM_VERIFY(ismatis,
"PetscBDDCSolver needs the matrix in unassembled format");
2791 ParFiniteElementSpace *fespace = opts.fespace;
2792 if (opts.netflux && !fespace)
2794 MFEM_WARNING(
"Don't know how to compute an auxiliary quadrature form without a ParFiniteElementSpace");
2801 ierr = MatISGetLocalMat(pA,&lA); CCHKERRQ(comm,ierr);
2802 ierr = MatNullSpaceCreate(PetscObjectComm((
PetscObject)lA),PETSC_TRUE,0,NULL,
2803 &nnsp); CCHKERRQ(PETSC_COMM_SELF,ierr);
2804 ierr = MatSetNearNullSpace(lA,nnsp); CCHKERRQ(PETSC_COMM_SELF,ierr);
2805 ierr = MatNullSpaceDestroy(&nnsp); CCHKERRQ(PETSC_COMM_SELF,ierr);
2809 ierr = PCSetType(pc,PCBDDC); PCHKERRQ(
obj,ierr);
2812 ierr = MatGetOwnershipRange(pA,&rst,NULL); PCHKERRQ(pA,ierr);
2813 ierr = MatGetLocalSize(pA,&nl,NULL); PCHKERRQ(pA,ierr);
2822 int vdim = fespace->GetVDim();
2829 ierr = MatSetBlockSize(pA,vdim); PCHKERRQ(pA,ierr);
2830 ierr = MatISGetLocalMat(pA,&lA); CCHKERRQ(PETSC_COMM_SELF,ierr);
2831 ierr = MatSetBlockSize(lA,vdim); PCHKERRQ(pA,ierr);
2840 ierr = PetscMalloc1(nf,&fields); CCHKERRQ(PETSC_COMM_SELF,ierr);
2853 ierr = ISCreateStride(comm,nlf,st,bs,&fields[i]); CCHKERRQ(comm,ierr);
2859 const FiniteElementCollection *fec = fespace->FEColl();
2860 bool h1space =
dynamic_cast<const H1_FECollection*
>(fec);
2863 ParFiniteElementSpace *fespace_coords = fespace;
2865 sdim = fespace->GetParMesh()->SpaceDimension();
2868 fespace_coords =
new ParFiniteElementSpace(fespace->GetParMesh(),fec,sdim,
2871 VectorFunctionCoefficient coeff_coords(sdim, func_coords);
2872 ParGridFunction gf_coords(fespace_coords);
2873 gf_coords.ProjectCoefficient(coeff_coords);
2874 HypreParVector *hvec_coords = gf_coords.ParallelProject();
2881 Vec pvec_coords,lvec_coords;
2882 ISLocalToGlobalMapping l2g;
2888 ierr = VecCreateMPIWithArray(comm,sdim,hvec_coords->Size(),
2889 hvec_coords->GlobalSize(),hvec_coords->GetData(),&pvec_coords);
2890 CCHKERRQ(comm,ierr);
2891 ierr = MatISGetLocalMat(pA,&lA); CCHKERRQ(PETSC_COMM_SELF,ierr);
2892 ierr = MatCreateVecs(lA,&lvec_coords,NULL); CCHKERRQ(PETSC_COMM_SELF,ierr);
2893 ierr = VecSetBlockSize(lvec_coords,sdim); CCHKERRQ(PETSC_COMM_SELF,ierr);
2894 ierr = MatGetLocalToGlobalMapping(pA,&l2g,NULL); CCHKERRQ(comm,ierr);
2895 ierr = MatGetLayouts(pA,&rmap,NULL); CCHKERRQ(comm,ierr);
2896 ierr = PetscSFCreate(comm,&sf); CCHKERRQ(comm,ierr);
2897 ierr = ISLocalToGlobalMappingGetIndices(l2g,&gidxs); CCHKERRQ(comm,ierr);
2898 ierr = ISLocalToGlobalMappingGetSize(l2g,&nleaves); CCHKERRQ(comm,ierr);
2899 ierr = PetscSFSetGraphLayout(sf,rmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs);
2900 CCHKERRQ(comm,ierr);
2901 ierr = ISLocalToGlobalMappingRestoreIndices(l2g,&gidxs); CCHKERRQ(comm,ierr);
2906 ierr = VecGetArrayRead(pvec_coords,&garray); CCHKERRQ(PETSC_COMM_SELF,ierr);
2907 ierr = VecGetArray(lvec_coords,&larray); CCHKERRQ(PETSC_COMM_SELF,ierr);
2908 ierr = PetscSFBcastBegin(sf,MPIU_SCALAR,garray,larray); CCHKERRQ(comm,ierr);
2909 ierr = PetscSFBcastEnd(sf,MPIU_SCALAR,garray,larray); CCHKERRQ(comm,ierr);
2910 ierr = VecRestoreArrayRead(pvec_coords,&garray); CCHKERRQ(PETSC_COMM_SELF,ierr);
2911 ierr = VecRestoreArray(lvec_coords,&larray); CCHKERRQ(PETSC_COMM_SELF,ierr);
2913 ierr = VecDestroy(&pvec_coords); CCHKERRQ(comm,ierr);
2914 ierr = MatNullSpaceCreateRigidBody(lvec_coords,&nnsp);
2915 CCHKERRQ(PETSC_COMM_SELF,ierr);
2916 ierr = VecDestroy(&lvec_coords); CCHKERRQ(PETSC_COMM_SELF,ierr);
2917 ierr = MatSetNearNullSpace(lA,nnsp); CCHKERRQ(PETSC_COMM_SELF,ierr);
2918 ierr = MatNullSpaceDestroy(&nnsp); CCHKERRQ(PETSC_COMM_SELF,ierr);
2919 ierr = PetscSFDestroy(&sf); CCHKERRQ(PETSC_COMM_SELF,ierr);
2923 ierr = PetscMalloc1(nl*sdim,&coords); CCHKERRQ(PETSC_COMM_SELF,ierr);
2933 ierr = ISGetLocalSize(fields[i],&nn); CCHKERRQ(comm,ierr);
2934 ierr = ISGetIndices(fields[i],&idxs); CCHKERRQ(comm,ierr);
2938 for (
PetscInt d = 0; d < sdim; d++)
2940 coords[sdim*idx+d] = data[sdim*j+d];
2943 ierr = ISRestoreIndices(fields[i],&idxs); CCHKERRQ(comm,ierr);
2948 ierr = PetscMemcpy(coords,hvec_coords->GetData(),nl*sdim*
sizeof(
PetscReal));
2949 CCHKERRQ(comm,ierr);
2951 if (fespace_coords != fespace)
2953 delete fespace_coords;
2960 IS dir = NULL, neu = NULL;
2963 Array<Mat> *l2l = NULL;
2964 if (opts.ess_dof_local || opts.nat_dof_local)
2969 MFEM_VERIFY(c,
"Local-to-local PETSc container not present");
2970 ierr = PetscContainerGetPointer(c,(
void**)&l2l); PCHKERRQ(c,ierr);
2977 PetscBool lpr = PETSC_FALSE,pr;
2978 if (opts.ess_dof) { lpr = PETSC_TRUE; }
2979 ierr = MPI_Allreduce(&lpr,&pr,1,MPIU_BOOL,MPI_LOR,comm);
2981 MFEM_VERIFY(lpr == pr,
"ess_dof should be collectively set");
2983 if (opts.nat_dof) { lpr = PETSC_TRUE; }
2984 ierr = MPI_Allreduce(&lpr,&pr,1,MPIU_BOOL,MPI_LOR,comm);
2986 MFEM_VERIFY(lpr == pr,
"nat_dof should be collectively set");
2989 ms[0] = -nf; ms[1] = nf;
2990 ierr = MPI_Allreduce(&ms,&Ms,2,MPIU_INT,MPI_MAX,comm);
2992 MFEM_VERIFY(-Ms[0] == Ms[1],
2993 "number of fields should be the same across processes");
3000 PetscInt st = opts.ess_dof_local ? 0 : rst;
3001 if (!opts.ess_dof_local)
3004 ierr = Convert_Array_IS(comm,
true,opts.ess_dof,st,&dir);
3005 CCHKERRQ(comm,ierr);
3006 ierr = PCBDDCSetDirichletBoundaries(pc,dir); PCHKERRQ(pc,ierr);
3011 ierr = Convert_Vmarks_IS(comm,*l2l,opts.ess_dof,st,&dir);
3012 CCHKERRQ(comm,ierr);
3013 ierr = PCBDDCSetDirichletBoundariesLocal(pc,dir); PCHKERRQ(pc,ierr);
3018 PetscInt st = opts.nat_dof_local ? 0 : rst;
3019 if (!opts.nat_dof_local)
3022 ierr = Convert_Array_IS(comm,
true,opts.nat_dof,st,&neu);
3023 CCHKERRQ(comm,ierr);
3024 ierr = PCBDDCSetNeumannBoundaries(pc,neu); PCHKERRQ(pc,ierr);
3029 ierr = Convert_Vmarks_IS(comm,*l2l,opts.nat_dof,st,&neu);
3030 CCHKERRQ(comm,ierr);
3031 ierr = PCBDDCSetNeumannBoundariesLocal(pc,neu); PCHKERRQ(pc,ierr);
3038 ierr = PCBDDCSetDofsSplitting(pc,nf,fields); PCHKERRQ(pc,ierr);
3042 ierr = ISDestroy(&fields[i]); CCHKERRQ(comm,ierr);
3044 ierr = PetscFree(fields); CCHKERRQ(PETSC_COMM_SELF,ierr);
3049 ierr = PCSetCoordinates(pc,sdim,nl,coords); PCHKERRQ(pc,ierr);
3051 ierr = PetscFree(coords); CCHKERRQ(PETSC_COMM_SELF,ierr);
3060 const FiniteElementCollection *fec = fespace->FEColl();
3061 bool edgespace, rtspace, h1space;
3062 bool needint = opts.netflux;
3063 bool tracespace, rt_tracespace, edge_tracespace;
3065 PetscBool B_is_Trans = PETSC_FALSE;
3067 ParMesh *pmesh = (ParMesh *) fespace->GetMesh();
3068 dim = pmesh->Dimension();
3069 vdim = fespace->GetVDim();
3070 h1space =
dynamic_cast<const H1_FECollection*
>(fec);
3071 rtspace =
dynamic_cast<const RT_FECollection*
>(fec);
3072 edgespace =
dynamic_cast<const ND_FECollection*
>(fec);
3073 edge_tracespace =
dynamic_cast<const ND_Trace_FECollection*
>(fec);
3074 rt_tracespace =
dynamic_cast<const RT_Trace_FECollection*
>(fec);
3075 tracespace = edge_tracespace || rt_tracespace;
3078 if (fespace->GetNE() > 0)
3082 p = fespace->GetOrder(0);
3086 p = fespace->GetFaceOrder(0);
3087 if (dim == 2) { p++; }
3098 MFEM_WARNING(
"Tracespace case doesn't work for H(curl) and p=2,"
3099 " not using auxiliary quadrature");
3105 FiniteElementCollection *vfec;
3108 vfec =
new H1_Trace_FECollection(p,dim);
3112 vfec =
new H1_FECollection(p,dim);
3114 ParFiniteElementSpace *vfespace =
new ParFiniteElementSpace(pmesh,vfec);
3115 ParDiscreteLinearOperator *grad;
3116 grad =
new ParDiscreteLinearOperator(vfespace,fespace);
3119 grad->AddTraceFaceInterpolator(
new GradientInterpolator);
3123 grad->AddDomainInterpolator(
new GradientInterpolator);
3127 HypreParMatrix *hG = grad->ParallelAssemble();
3128 PetscParMatrix *G =
new PetscParMatrix(hG,
PETSC_MATAIJ);
3132 PetscBool conforming = PETSC_TRUE;
3133 if (pmesh->Nonconforming()) { conforming = PETSC_FALSE; }
3134 ierr = PCBDDCSetDiscreteGradient(pc,*G,p,0,PETSC_TRUE,conforming);
3146 MFEM_WARNING(
"Tracespace case doesn't work for H(div), not using"
3147 " auxiliary quadrature");
3153 if (vdim != dim) { needint =
false; }
3156 PetscParMatrix *
B = NULL;
3162 FiniteElementCollection *auxcoll;
3163 if (tracespace) { auxcoll =
new RT_Trace_FECollection(p,dim); }
3168 auxcoll =
new H1_FECollection(std::max(p-1,1),dim);
3172 auxcoll =
new L2_FECollection(p,dim);
3175 ParFiniteElementSpace *pspace =
new ParFiniteElementSpace(pmesh,auxcoll);
3176 ParMixedBilinearForm *b =
new ParMixedBilinearForm(fespace,pspace);
3182 b->AddTraceFaceIntegrator(
new VectorFECurlIntegrator);
3186 b->AddDomainIntegrator(
new VectorFECurlIntegrator);
3193 b->AddTraceFaceIntegrator(
new VectorFEDivergenceIntegrator);
3197 b->AddDomainIntegrator(
new VectorFEDivergenceIntegrator);
3202 b->AddDomainIntegrator(
new VectorDivergenceIntegrator);
3207 b->ParallelAssemble(Bh);
3209 Bh.SetOperatorOwner(
false);
3214 ierr = MatTranspose(pB,MAT_INPLACE_MATRIX,&pB); PCHKERRQ(pA,ierr);
3215 if (!opts.ess_dof_local)
3217 ierr = MatZeroRowsIS(pB,dir,0.,NULL,NULL); PCHKERRQ(pA,ierr);
3221 ierr = MatZeroRowsLocalIS(pB,dir,0.,NULL,NULL); PCHKERRQ(pA,ierr);
3223 B_is_Trans = PETSC_TRUE;
3232 ierr = PCBDDCSetDivergenceMat(pc,*B,B_is_Trans,NULL); PCHKERRQ(pc,ierr);
3236 ierr = ISDestroy(&dir); PCHKERRQ(pc,ierr);
3237 ierr = ISDestroy(&neu); PCHKERRQ(pc,ierr);
3242 const std::string &prefix)
3245 BDDCSolverConstructor(opts);
3251 const std::string &prefix)
3254 BDDCSolverConstructor(opts);
3259 const string &prefix)
3263 ierr = PCSetType(pc,PCFIELDSPLIT); PCHKERRQ(pc,ierr);
3266 ierr = PCGetOperators(pc,&pA,NULL); PCHKERRQ(pc,ierr);
3270 ierr = PetscObjectTypeCompare((
PetscObject)pA,MATNEST,&isnest);
3276 ierr = MatNestGetSize(pA,&nr,NULL); PCHKERRQ(pc,ierr);
3277 ierr = PetscCalloc1(nr,&isrow); CCHKERRQ(PETSC_COMM_SELF,ierr);
3278 ierr = MatNestGetISs(pA,isrow,NULL); PCHKERRQ(pc,ierr);
3288 ierr = PCFieldSplitSetIS(pc,NULL,isrow[i]); PCHKERRQ(pc,ierr);
3290 ierr = PetscFree(isrow); CCHKERRQ(PETSC_COMM_SELF,ierr);
3296 const std::string &prefix)
3301 ierr = SNESCreate(comm, &snes); CCHKERRQ(comm, ierr);
3303 ierr = PetscObjectGetClassId(
obj, &
cid); PCHKERRQ(
obj, ierr);
3304 ierr = SNESSetOptionsPrefix(snes, prefix.c_str()); PCHKERRQ(snes, ierr);
3311 const std::string &prefix)
3316 ierr = SNESCreate(comm, &snes); CCHKERRQ(comm, ierr);
3318 ierr = PetscObjectGetClassId(
obj, &
cid); PCHKERRQ(
obj, ierr);
3319 ierr = SNESSetOptionsPrefix(snes, prefix.c_str()); PCHKERRQ(snes, ierr);
3331 ierr = PetscObjectGetComm(
obj,&comm); PCHKERRQ(
obj, ierr);
3332 ierr = SNESDestroy(&snes); CCHKERRQ(comm, ierr);
3344 ierr = SNESGetFunction(snes, NULL, NULL, &fctx);
3345 PCHKERRQ(snes, ierr);
3346 ierr = SNESGetJacobian(snes, NULL, NULL, NULL, &jctx);
3347 PCHKERRQ(snes, ierr);
3350 (
void*)&op == fctx &&
3351 (
void*)&op == jctx);
3352 ierr = MPI_Allreduce(&ls,&gs,1,MPIU_BOOL,MPI_LAND,
3354 PCHKERRQ(snes,ierr);
3357 ierr = SNESReset(snes); PCHKERRQ(snes,ierr);
3368 ierr = SNESGetLineSearch(snes, &ls); PCHKERRQ(snes,ierr);
3369 ierr = SNESLineSearchSetType(ls, SNESLINESEARCHBT); PCHKERRQ(snes,ierr);
3372 __mfem_snes_ctx *snes_ctx = (__mfem_snes_ctx*)
private_ctx;
3374 ierr = SNESSetFunction(snes, NULL, __mfem_snes_function, (
void *)snes_ctx);
3375 PCHKERRQ(snes, ierr);
3376 ierr = SNESSetJacobian(snes, NULL, NULL, __mfem_snes_jacobian,
3378 PCHKERRQ(snes, ierr);
3390 __mfem_snes_ctx *snes_ctx = (__mfem_snes_ctx*)
private_ctx;
3391 snes_ctx->jacType = jacType;
3397 __mfem_snes_ctx *snes_ctx = (__mfem_snes_ctx*)
private_ctx;
3398 snes_ctx->objective = objfn;
3401 ierr = SNESSetObjective(snes, __mfem_snes_objective, (
void *)snes_ctx);
3402 PCHKERRQ(snes, ierr);
3409 __mfem_snes_ctx *snes_ctx = (__mfem_snes_ctx*)
private_ctx;
3410 snes_ctx->postcheck = post;
3414 ierr = SNESGetLineSearch(snes, &ls); PCHKERRQ(snes,ierr);
3415 ierr = SNESLineSearchSetPostCheck(ls, __mfem_snes_postcheck, (
void *)snes_ctx);
3425 __mfem_snes_ctx *snes_ctx = (__mfem_snes_ctx*)
private_ctx;
3426 snes_ctx->update = update;
3429 ierr = SNESSetUpdate(snes, __mfem_snes_update); PCHKERRQ(snes, ierr);
3436 bool b_nonempty = b.
Size();
3440 if (b_nonempty) { B->PlaceArray(b.
GetData()); }
3448 ierr = SNESSolve(snes, B->x, X->x); PCHKERRQ(snes, ierr);
3450 if (b_nonempty) { B->ResetArray(); }
3460 ierr = TSCreate(comm,&ts); CCHKERRQ(comm,ierr);
3462 ierr = PetscObjectGetClassId(
obj,&
cid); PCHKERRQ(
obj,ierr);
3463 ierr = TSSetOptionsPrefix(ts, prefix.c_str()); PCHKERRQ(ts, ierr);
3469 ierr = TSSetMaxSteps(ts,PETSC_MAX_INT-1);
3471 ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
3474 ierr = TSGetAdapt(ts,&tsad);
3476 ierr = TSAdaptSetType(tsad,TSADAPTNONE);
3484 ierr = PetscObjectGetComm(
obj,&comm); PCHKERRQ(
obj,ierr);
3485 ierr = TSDestroy(&ts); CCHKERRQ(comm,ierr);
3493 __mfem_ts_ctx *ts_ctx = (__mfem_ts_ctx*)
private_ctx;
3496 ierr = TSReset(ts); PCHKERRQ(ts,ierr);
3499 ts_ctx->cached_shift = std::numeric_limits<PetscReal>::min();
3500 ts_ctx->cached_ijacstate = -1;
3501 ts_ctx->cached_rhsjacstate = -1;
3502 ts_ctx->cached_splits_xstate = -1;
3503 ts_ctx->cached_splits_xdotstate = -1;
3511 ierr = TSSetIFunction(ts, NULL, __mfem_ts_ifunction, (
void *)ts_ctx);
3513 ierr = TSSetIJacobian(ts, NULL, NULL, __mfem_ts_ijacobian, (
void *)ts_ctx);
3515 ierr = TSSetEquationType(ts, TS_EQ_IMPLICIT);
3522 ierr = TSSetEquationType(ts, TS_EQ_EXPLICIT);
3525 ierr = TSSetRHSFunction(ts, NULL, __mfem_ts_rhsfunction, (
void *)ts_ctx);
3527 ierr = TSSetRHSJacobian(ts, NULL, NULL, __mfem_ts_rhsjacobian, (
void *)ts_ctx);
3536 ierr = TSSetSolution(ts,X); PCHKERRQ(ts,ierr);
3539 PetscBool use = PETSC_TRUE;
3540 ierr = PetscOptionsGetBool(NULL,NULL,
"-mfem_use_splitjac",&use,NULL);
3543 ierr = PetscObjectComposeFunction((
PetscObject)ts,
"TSComputeSplitJacobians_C",
3544 __mfem_ts_computesplits);
3549 ierr = PetscObjectComposeFunction((
PetscObject)ts,
"TSComputeSplitJacobians_C",
3557 __mfem_ts_ctx *ts_ctx = (__mfem_ts_ctx*)
private_ctx;
3558 ts_ctx->jacType = jacType;
3563 __mfem_ts_ctx *ts_ctx = (__mfem_ts_ctx*)
private_ctx;
3564 return ts_ctx->type;
3569 __mfem_ts_ctx *ts_ctx = (__mfem_ts_ctx*)
private_ctx;
3572 ts_ctx->type = type;
3575 ierr = TSSetProblemType(ts, TS_LINEAR);
3580 ierr = TSSetProblemType(ts, TS_NONLINEAR);
3589 ierr = TSSetTime(ts, t); PCHKERRQ(ts, ierr);
3590 ierr = TSSetTimeStep(ts, dt); PCHKERRQ(ts, ierr);
3598 ierr = TSSetSolution(ts, *X); PCHKERRQ(ts, ierr);
3599 ierr = TSStep(ts); PCHKERRQ(ts, ierr);
3605 ierr = TSGetTime(ts,&pt); PCHKERRQ(ts,ierr);
3614 ierr = TSSetTime(ts, t); PCHKERRQ(ts, ierr);
3615 ierr = TSSetTimeStep(ts, dt); PCHKERRQ(ts, ierr);
3616 ierr = TSSetMaxTime(ts, t_final); PCHKERRQ(ts, ierr);
3617 ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);
3629 __mfem_ts_ctx *ts_ctx = (__mfem_ts_ctx*)
private_ctx;
3630 ts_ctx->cached_shift = std::numeric_limits<PetscReal>::min();
3631 ts_ctx->cached_ijacstate = -1;
3632 ts_ctx->cached_rhsjacstate = -1;
3633 ts_ctx->cached_splits_xstate = -1;
3634 ts_ctx->cached_splits_xdotstate = -1;
3637 ierr = TSSolve(ts, X->x); PCHKERRQ(ts, ierr);
3642 ierr = TSGetTime(ts, &pt); PCHKERRQ(ts,ierr);
3644 ierr = TSGetTimeStep(ts,&pt); PCHKERRQ(ts,ierr);
3650 #include "petsc/private/petscimpl.h"
3651 #include "petsc/private/matimpl.h"
3657 __mfem_monitor_ctx *monctx = (__mfem_monitor_ctx*)ctx;
3659 PetscFunctionBeginUser;
3662 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,
"Missing monitor context");
3674 PetscFunctionReturn(0);
3677 static PetscErrorCode __mfem_ts_ifunction(
TS ts,
PetscReal t, Vec x, Vec xp,
3680 __mfem_ts_ctx* ts_ctx = (__mfem_ts_ctx*)ctx;
3681 PetscErrorCode ierr;
3683 PetscFunctionBeginUser;
3691 if (ts_ctx->bchandler)
3696 if (!ts_ctx->work) { ts_ctx->work =
new mfem::Vector(xx.Size()); }
3697 if (!ts_ctx->work2) { ts_ctx->work2 =
new mfem::Vector(xx.Size()); }
3703 bchandler->
ZeroBC(yy,*txp);
3715 ierr = PetscObjectStateIncrease((
PetscObject)f); CHKERRQ(ierr);
3716 PetscFunctionReturn(0);
3719 static PetscErrorCode __mfem_ts_rhsfunction(
TS ts,
PetscReal t, Vec x, Vec f,
3722 __mfem_ts_ctx* ts_ctx = (__mfem_ts_ctx*)ctx;
3723 PetscErrorCode ierr;
3725 PetscFunctionBeginUser;
3726 if (ts_ctx->bchandler) { MFEM_ABORT(
"RHS evaluation with bc not implemented"); }
3736 ierr = PetscObjectStateIncrease((
PetscObject)f); CHKERRQ(ierr);
3737 PetscFunctionReturn(0);
3740 static PetscErrorCode __mfem_ts_ijacobian(
TS ts,
PetscReal t, Vec x,
3744 __mfem_ts_ctx* ts_ctx = (__mfem_ts_ctx*)ctx;
3749 PetscObjectState state;
3750 PetscErrorCode ierr;
3752 PetscFunctionBeginUser;
3756 ierr = PetscObjectStateGet((
PetscObject)A,&state); CHKERRQ(ierr);
3758 std::abs(ts_ctx->cached_shift/shift - 1.0) < eps &&
3759 state == ts_ctx->cached_ijacstate) { PetscFunctionReturn(0); }
3766 ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr);
3767 ierr = VecGetArrayRead(xp,(
const PetscScalar**)&array); CHKERRQ(ierr);
3769 ierr = VecRestoreArrayRead(xp,(
const PetscScalar**)&array); CHKERRQ(ierr);
3770 ierr = VecGetArrayRead(x,(
const PetscScalar**)&array); CHKERRQ(ierr);
3771 if (!ts_ctx->bchandler)
3778 if (!ts_ctx->work) { ts_ctx->work =
new mfem::Vector(n); }
3785 ierr = VecRestoreArrayRead(x,(
const PetscScalar**)&array); CHKERRQ(ierr);
3789 if (!ts_ctx->bchandler) {
delete xx; }
3790 ts_ctx->cached_shift = shift;
3793 bool delete_pA =
false;
3797 pA->
GetType() != ts_ctx->jacType))
3805 if (ts_ctx->bchandler)
3813 PetscObjectState nonzerostate;
3814 ierr = MatGetNonzeroState(P,&nonzerostate); CHKERRQ(ierr);
3819 ierr = MatHeaderReplace(P,&B); CHKERRQ(ierr);
3820 if (delete_pA) {
delete pA; }
3825 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
3826 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
3839 ierr = PetscObjectTypeCompare((
PetscObject)P,MATNEST,&isnest);
3841 if (isnest) { P->nonzerostate = nonzerostate + 1; }
3844 ierr = PetscObjectStateGet((
PetscObject)P,&ts_ctx->cached_ijacstate);
3850 ierr = MatGetType(P,&mtype); CHKERRQ(ierr);
3851 ierr = TSGetDM(ts,&dm); CHKERRQ(ierr);
3852 ierr = DMSetMatType(dm,mtype); CHKERRQ(ierr);
3853 ierr = DMShellSetMatrix(dm,P); CHKERRQ(ierr);
3854 PetscFunctionReturn(0);
3857 static PetscErrorCode __mfem_ts_computesplits(
TS ts,
PetscReal t,Vec x,Vec xp,
3861 __mfem_ts_ctx* ts_ctx;
3865 PetscObjectState state;
3866 PetscBool rx = PETSC_TRUE, rxp = PETSC_TRUE;
3867 PetscBool assembled;
3868 PetscErrorCode ierr;
3870 PetscFunctionBeginUser;
3871 ierr = TSGetIJacobian(ts,NULL,NULL,NULL,(
void**)&ts_ctx); CHKERRQ(ierr);
3874 ierr = PetscObjectStateGet((
PetscObject)Jx,&state); CHKERRQ(ierr);
3876 state == ts_ctx->cached_splits_xstate) { rx = PETSC_FALSE; }
3877 ierr = PetscObjectStateGet((
PetscObject)Jxp,&state); CHKERRQ(ierr);
3879 state == ts_ctx->cached_splits_xdotstate) { rxp = PETSC_FALSE; }
3880 if (!rx && !rxp) { PetscFunctionReturn(0); }
3887 ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr);
3888 ierr = VecGetArrayRead(xp,(
const PetscScalar**)&array); CHKERRQ(ierr);
3890 ierr = VecRestoreArrayRead(xp,(
const PetscScalar**)&array); CHKERRQ(ierr);
3891 ierr = VecGetArrayRead(x,(
const PetscScalar**)&array); CHKERRQ(ierr);
3892 if (!ts_ctx->bchandler)
3899 if (!ts_ctx->work) { ts_ctx->work =
new mfem::Vector(n); }
3906 ierr = VecRestoreArrayRead(x,(
const PetscScalar**)&array); CHKERRQ(ierr);
3915 bool delete_mat =
false;
3919 pJx->
GetType() != ts_ctx->jacType))
3924 ierr = PetscObjectReference((
PetscObject)B); CHKERRQ(ierr);
3932 ierr = MatAssembled(Jx,&assembled); CHKERRQ(ierr);
3935 ierr = MatCopy(*pJx,Jx,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
3940 ierr = MatDuplicate(*pJx,MAT_COPY_VALUES,&B); CHKERRQ(ierr);
3941 ierr = MatHeaderReplace(Jx,&B); CHKERRQ(ierr);
3944 if (delete_mat) {
delete pJx; }
3948 if (ts_ctx->bchandler)
3965 pJxp->
GetType() != ts_ctx->jacType))
3970 ierr = PetscObjectReference((
PetscObject)B); CHKERRQ(ierr);
3973 &oJxp,ts_ctx->jacType);
3977 ierr = MatAssembled(Jxp,&assembled); CHKERRQ(ierr);
3980 ierr = MatCopy(*pJxp,Jxp,SAME_NONZERO_PATTERN); CHKERRQ(ierr);
3985 ierr = MatDuplicate(*pJxp,MAT_COPY_VALUES,&B); CHKERRQ(ierr);
3986 ierr = MatHeaderReplace(Jxp,&B); CHKERRQ(ierr);
3988 if (delete_mat) {
delete pJxp; }
3992 if (ts_ctx->bchandler)
4003 ierr = MatAXPY(*pJxp,-1.0,*pJx,SAME_NONZERO_PATTERN); PCHKERRQ(ts,ierr);
4009 ierr = MatAssemblyBegin(Ax,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
4010 ierr = MatAssemblyEnd(Ax,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
4012 if (Axp && Axp != Jxp)
4014 ierr = MatAssemblyBegin(Axp,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
4015 ierr = MatAssemblyEnd(Axp,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
4019 ierr = PetscObjectStateGet((
PetscObject)Jx,&ts_ctx->cached_splits_xstate);
4021 ierr = PetscObjectStateGet((
PetscObject)Jxp,&ts_ctx->cached_splits_xdotstate);
4026 if (!ts_ctx->bchandler) {
delete xx; }
4027 PetscFunctionReturn(0);
4030 static PetscErrorCode __mfem_ts_rhsjacobian(
TS ts,
PetscReal t, Vec x,
4033 __mfem_ts_ctx* ts_ctx = (__mfem_ts_ctx*)ctx;
4037 PetscObjectState state;
4038 PetscErrorCode ierr;
4040 PetscFunctionBeginUser;
4042 ierr = PetscObjectStateGet((
PetscObject)A,&state); CHKERRQ(ierr);
4044 state == ts_ctx->cached_rhsjacstate) { PetscFunctionReturn(0); }
4051 ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr);
4052 ierr = VecGetArrayRead(x,(
const PetscScalar**)&array); CHKERRQ(ierr);
4053 if (!ts_ctx->bchandler)
4060 if (!ts_ctx->work) { ts_ctx->work =
new mfem::Vector(n); }
4067 ierr = VecRestoreArrayRead(x,(
const PetscScalar**)&array); CHKERRQ(ierr);
4071 if (!ts_ctx->bchandler) {
delete xx; }
4074 bool delete_pA =
false;
4078 pA->
GetType() != ts_ctx->jacType))
4086 if (ts_ctx->bchandler)
4094 PetscObjectState nonzerostate;
4095 ierr = MatGetNonzeroState(P,&nonzerostate); CHKERRQ(ierr);
4100 ierr = MatHeaderReplace(P,&B); CHKERRQ(ierr);
4101 if (delete_pA) {
delete pA; }
4113 ierr = PetscObjectTypeCompare((
PetscObject)P,MATNEST,&isnest);
4115 if (isnest) { P->nonzerostate = nonzerostate + 1; }
4120 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
4121 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
4127 ierr = TSRHSJacobianSetReuse(ts,PETSC_TRUE); PCHKERRQ(ts,ierr);
4129 ierr = PetscObjectStateGet((
PetscObject)A,&ts_ctx->cached_rhsjacstate);
4135 ierr = MatGetType(P,&mtype); CHKERRQ(ierr);
4136 ierr = TSGetDM(ts,&dm); CHKERRQ(ierr);
4137 ierr = DMSetMatType(dm,mtype); CHKERRQ(ierr);
4138 ierr = DMShellSetMatrix(dm,P); CHKERRQ(ierr);
4139 PetscFunctionReturn(0);
4145 __mfem_monitor_ctx *monctx = (__mfem_monitor_ctx*)ctx;
4147 PetscFunctionBeginUser;
4150 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,
"Missing monitor context");
4159 PetscErrorCode ierr;
4161 ierr = SNESGetSolution(snes,&x); CHKERRQ(ierr);
4168 PetscErrorCode ierr;
4170 ierr = SNESGetFunction(snes,&x,NULL,NULL); CHKERRQ(ierr);
4175 PetscFunctionReturn(0);
4178 static PetscErrorCode __mfem_snes_jacobian(
SNES snes, Vec x,
Mat A,
Mat P,
4183 PetscErrorCode ierr;
4185 __mfem_snes_ctx *snes_ctx = (__mfem_snes_ctx*)ctx;
4187 PetscFunctionBeginUser;
4188 ierr = VecGetArrayRead(x,(
const PetscScalar**)&array); CHKERRQ(ierr);
4189 ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr);
4190 if (!snes_ctx->bchandler)
4197 if (!snes_ctx->work) { snes_ctx->work =
new mfem::Vector(n); }
4200 xx = snes_ctx->work;
4206 ierr = VecRestoreArrayRead(x,(
const PetscScalar**)&array); CHKERRQ(ierr);
4207 if (!snes_ctx->bchandler) {
delete xx; }
4210 bool delete_pA =
false;
4214 pA->
GetType() != snes_ctx->jacType))
4222 if (snes_ctx->bchandler)
4230 PetscObjectState nonzerostate;
4231 ierr = MatGetNonzeroState(P,&nonzerostate); CHKERRQ(ierr);
4235 ierr = MatHeaderReplace(P,&B); CHKERRQ(ierr);
4236 if (delete_pA) {
delete pA; }
4248 ierr = PetscObjectTypeCompare((
PetscObject)P,MATNEST,&isnest);
4250 if (isnest) { P->nonzerostate = nonzerostate + 1; }
4255 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
4256 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
4262 ierr = MatGetType(P,&mtype); CHKERRQ(ierr);
4263 ierr = SNESGetDM(snes,&dm); CHKERRQ(ierr);
4264 ierr = DMSetMatType(dm,mtype); CHKERRQ(ierr);
4265 ierr = DMShellSetMatrix(dm,P); CHKERRQ(ierr);
4266 PetscFunctionReturn(0);
4269 static PetscErrorCode __mfem_snes_function(
SNES snes, Vec x, Vec f,
void *ctx)
4271 __mfem_snes_ctx* snes_ctx = (__mfem_snes_ctx*)ctx;
4273 PetscFunctionBeginUser;
4276 if (snes_ctx->bchandler)
4283 snes_ctx->op->Mult(*txx,ff);
4290 snes_ctx->op->Mult(xx,ff);
4293 ierr = PetscObjectStateIncrease((
PetscObject)f); CHKERRQ(ierr);
4294 PetscFunctionReturn(0);
4297 static PetscErrorCode __mfem_snes_objective(
SNES snes, Vec x,
PetscReal *f,
4300 __mfem_snes_ctx* snes_ctx = (__mfem_snes_ctx*)ctx;
4302 PetscFunctionBeginUser;
4303 if (!snes_ctx->objective)
4305 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,
"Missing objective function");
4309 (*snes_ctx->objective)(snes_ctx->op,xx,&lf);
4311 PetscFunctionReturn(0);
4314 static PetscErrorCode __mfem_snes_postcheck(SNESLineSearch ls,Vec X,Vec Y,Vec W,
4315 PetscBool *cy,PetscBool *cw,
void* ctx)
4317 __mfem_snes_ctx* snes_ctx = (__mfem_snes_ctx*)ctx;
4318 bool lcy =
false,lcw =
false;
4320 PetscFunctionBeginUser;
4324 (*snes_ctx->postcheck)(snes_ctx->op,x,y,w,lcy,lcw);
4325 if (lcy) { *cy = PETSC_TRUE; }
4326 if (lcw) { *cw = PETSC_TRUE; }
4327 PetscFunctionReturn(0);
4330 static PetscErrorCode __mfem_snes_update(
SNES snes,
PetscInt it)
4333 __mfem_snes_ctx* snes_ctx;
4335 PetscFunctionBeginUser;
4337 ierr = SNESGetFunction(snes,&F,NULL,(
void **)&snes_ctx); CHKERRQ(ierr);
4338 ierr = SNESGetSolution(snes,&X); CHKERRQ(ierr);
4341 ierr = VecDuplicate(X,&pX); CHKERRQ(ierr);
4344 ierr = VecDestroy(&pX); CHKERRQ(ierr);
4348 if (!pX) SETERRQ(PetscObjectComm((
PetscObject)snes),PETSC_ERR_USER,
4349 "Missing previous solution");
4350 ierr = SNESGetSolutionUpdate(snes,&dX); CHKERRQ(ierr);
4355 (*snes_ctx->update)(snes_ctx->op,it,f,x,dx,px);
4357 ierr = VecCopy(X,pX); CHKERRQ(ierr);
4358 PetscFunctionReturn(0);
4364 __mfem_monitor_ctx *monctx = (__mfem_monitor_ctx*)ctx;
4366 PetscFunctionBeginUser;
4369 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,
"Missing monitor context");
4378 PetscErrorCode ierr;
4380 ierr = KSPBuildSolution(ksp,NULL,&x); CHKERRQ(ierr);
4387 PetscErrorCode ierr;
4389 ierr = KSPBuildResidual(ksp,NULL,NULL,&x); CHKERRQ(ierr);
4394 PetscFunctionReturn(0);
4397 static PetscErrorCode __mfem_mat_shell_apply(
Mat A, Vec x, Vec y)
4400 PetscErrorCode ierr;
4402 PetscFunctionBeginUser;
4403 ierr = MatShellGetContext(A,(
void **)&op); CHKERRQ(ierr);
4404 if (!op) { SETERRQ(PetscObjectComm((
PetscObject)A),PETSC_ERR_LIB,
"Missing operator"); }
4409 ierr = PetscObjectStateIncrease((
PetscObject)y); CHKERRQ(ierr);
4410 PetscFunctionReturn(0);
4413 static PetscErrorCode __mfem_mat_shell_apply_transpose(
Mat A, Vec x, Vec y)
4416 PetscErrorCode ierr;
4418 PetscFunctionBeginUser;
4419 ierr = MatShellGetContext(A,(
void **)&op); CHKERRQ(ierr);
4420 if (!op) { SETERRQ(PetscObjectComm((
PetscObject)A),PETSC_ERR_LIB,
"Missing operator"); }
4425 ierr = PetscObjectStateIncrease((
PetscObject)y); CHKERRQ(ierr);
4426 PetscFunctionReturn(0);
4429 static PetscErrorCode __mfem_mat_shell_copy(
Mat A,
Mat B, MatStructure str)
4432 PetscErrorCode ierr;
4434 PetscFunctionBeginUser;
4435 ierr = MatShellGetContext(A,(
void **)&op); CHKERRQ(ierr);
4436 if (!op) { SETERRQ(PetscObjectComm((
PetscObject)A),PETSC_ERR_LIB,
"Missing operator"); }
4437 ierr = MatShellSetContext(B,(
void *)op); CHKERRQ(ierr);
4438 PetscFunctionReturn(0);
4441 static PetscErrorCode __mfem_mat_shell_destroy(
Mat A)
4443 PetscFunctionBeginUser;
4444 PetscFunctionReturn(0);
4447 static PetscErrorCode __mfem_pc_shell_view(
PC pc, PetscViewer viewer)
4449 __mfem_pc_shell_ctx *ctx;
4450 PetscErrorCode ierr;
4452 PetscFunctionBeginUser;
4453 ierr = PCShellGetContext(pc,(
void **)&ctx); CHKERRQ(ierr);
4457 ierr = PetscObjectTypeCompare((
PetscObject)viewer,PETSCVIEWERASCII,&isascii);
4464 ierr = PCView(*ppc,viewer); CHKERRQ(ierr);
4470 ierr = PetscViewerASCIIPrintf(viewer,
4471 "No information available on the mfem::Solver\n");
4475 if (isascii && ctx->factory)
4477 ierr = PetscViewerASCIIPrintf(viewer,
4478 "Number of preconditioners created by the factory %lu\n",ctx->numprec);
4482 PetscFunctionReturn(0);
4485 static PetscErrorCode __mfem_pc_shell_apply(
PC pc, Vec x, Vec y)
4487 __mfem_pc_shell_ctx *ctx;
4488 PetscErrorCode ierr;
4490 PetscFunctionBeginUser;
4493 ierr = PCShellGetContext(pc,(
void **)&ctx); CHKERRQ(ierr);
4496 ctx->op->Mult(xx,yy);
4498 ierr = PetscObjectStateIncrease((
PetscObject)y); CHKERRQ(ierr);
4504 PetscFunctionReturn(0);
4507 static PetscErrorCode __mfem_pc_shell_apply_transpose(
PC pc, Vec x, Vec y)
4509 __mfem_pc_shell_ctx *ctx;
4510 PetscErrorCode ierr;
4512 PetscFunctionBeginUser;
4515 ierr = PCShellGetContext(pc,(
void **)&ctx); CHKERRQ(ierr);
4518 ctx->op->MultTranspose(xx,yy);
4520 ierr = PetscObjectStateIncrease((
PetscObject)y); CHKERRQ(ierr);
4526 PetscFunctionReturn(0);
4529 static PetscErrorCode __mfem_pc_shell_setup(
PC pc)
4531 __mfem_pc_shell_ctx *ctx;
4533 PetscFunctionBeginUser;
4534 ierr = PCShellGetContext(pc,(
void **)&ctx); CHKERRQ(ierr);
4545 ierr = PCGetOperators(pc,NULL,&B); CHKERRQ(ierr);
4554 PetscFunctionReturn(0);
4557 static PetscErrorCode __mfem_pc_shell_destroy(
PC pc)
4559 __mfem_pc_shell_ctx *ctx;
4560 PetscErrorCode ierr;
4562 PetscFunctionBeginUser;
4563 ierr = PCShellGetContext(pc,(
void **)&ctx); CHKERRQ(ierr);
4569 PetscFunctionReturn(0);
4572 static PetscErrorCode __mfem_array_container_destroy(
void *ptr)
4574 PetscErrorCode ierr;
4576 PetscFunctionBeginUser;
4577 ierr = PetscFree(ptr); CHKERRQ(ierr);
4578 PetscFunctionReturn(0);
4581 static PetscErrorCode __mfem_matarray_container_destroy(
void *ptr)
4584 PetscErrorCode ierr;
4586 PetscFunctionBeginUser;
4587 for (
int i=0; i<a->
Size(); i++)
4591 ierr = MatDestroy(&M); CCHKERRQ(comm,ierr);
4594 PetscFunctionReturn(0);
4597 static PetscErrorCode __mfem_monitor_ctx_destroy(
void **ctx)
4599 PetscErrorCode ierr;
4601 PetscFunctionBeginUser;
4602 ierr = PetscFree(*ctx); CHKERRQ(ierr);
4603 PetscFunctionReturn(0);
4608 PetscErrorCode MakeShellPC(
PC pc,
mfem::Solver &precond,
bool ownsop)
4610 PetscFunctionBeginUser;
4611 __mfem_pc_shell_ctx *ctx =
new __mfem_pc_shell_ctx;
4613 ctx->ownsop = ownsop;
4614 ctx->factory = NULL;
4620 ierr = PCSetType(pc,PCNONE); CHKERRQ(ierr);
4622 ierr = PCSetType(pc,PCSHELL); CHKERRQ(ierr);
4623 ierr = PCShellSetName(pc,
"MFEM Solver (unknown Pmat)"); CHKERRQ(ierr);
4624 ierr = PCShellSetContext(pc,(
void *)ctx); CHKERRQ(ierr);
4625 ierr = PCShellSetApply(pc,__mfem_pc_shell_apply); CHKERRQ(ierr);
4626 ierr = PCShellSetApplyTranspose(pc,__mfem_pc_shell_apply_transpose);
4628 ierr = PCShellSetSetUp(pc,__mfem_pc_shell_setup); CHKERRQ(ierr);
4629 ierr = PCShellSetView(pc,__mfem_pc_shell_view); CHKERRQ(ierr);
4630 ierr = PCShellSetDestroy(pc,__mfem_pc_shell_destroy); CHKERRQ(ierr);
4631 PetscFunctionReturn(0);
4636 PetscErrorCode MakeShellPCWithFactory(
PC pc,
4639 PetscFunctionBeginUser;
4640 __mfem_pc_shell_ctx *ctx =
new __mfem_pc_shell_ctx;
4643 ctx->factory = factory;
4649 ierr = PCSetType(pc,PCNONE); CHKERRQ(ierr);
4651 ierr = PCSetType(pc,PCSHELL); CHKERRQ(ierr);
4652 ierr = PCShellSetName(pc,factory->
GetName()); CHKERRQ(ierr);
4653 ierr = PCShellSetContext(pc,(
void *)ctx); CHKERRQ(ierr);
4654 ierr = PCShellSetApply(pc,__mfem_pc_shell_apply); CHKERRQ(ierr);
4655 ierr = PCShellSetApplyTranspose(pc,__mfem_pc_shell_apply_transpose);
4657 ierr = PCShellSetSetUp(pc,__mfem_pc_shell_setup); CHKERRQ(ierr);
4658 ierr = PCShellSetView(pc,__mfem_pc_shell_view); CHKERRQ(ierr);
4659 ierr = PCShellSetDestroy(pc,__mfem_pc_shell_destroy); CHKERRQ(ierr);
4660 PetscFunctionReturn(0);
4665 static PetscErrorCode Convert_Array_IS(MPI_Comm comm,
bool islist,
4670 const int *data = list ? list->
GetData() : NULL;
4671 PetscErrorCode ierr;
4673 PetscFunctionBeginUser;
4674 ierr = PetscMalloc1(n,&idxs); CHKERRQ(ierr);
4677 for (
PetscInt i=0; i<n; i++) { idxs[i] = data[i] + st; }
4684 if (data[i]) { idxs[cum++] = i+st; }
4688 ierr = ISCreateGeneral(comm,n,idxs,PETSC_OWN_POINTER,is);
4690 PetscFunctionReturn(0);
4696 static PetscErrorCode Convert_Vmarks_IS(MPI_Comm comm,
4704 PetscErrorCode ierr;
4706 PetscFunctionBeginUser;
4707 for (
int i = 0; i < pl2l.
Size(); i++)
4711 ierr = MatGetRowIJ(pl2l[i],0,PETSC_FALSE,PETSC_FALSE,&m,(
const PetscInt**)&ii,
4712 (
const PetscInt**)&jj,&done); CHKERRQ(ierr);
4713 MFEM_VERIFY(done,
"Unable to perform MatGetRowIJ on " << i <<
" l2l matrix");
4714 ierr = MatGetSize(pl2l[i],NULL,&n); CHKERRQ(ierr);
4715 #if defined(PETSC_USE_64BIT_INDICES)
4716 int nnz = (int)ii[m];
4717 int *mii =
new int[m+1];
4718 int *mjj =
new int[nnz];
4719 for (
int j = 0; j < m+1; j++) { mii[j] = (int)ii[j]; }
4720 for (
int j = 0; j < nnz; j++) { mjj[j] = (int)jj[j]; }
4725 ierr = MatRestoreRowIJ(pl2l[i],0,PETSC_FALSE,PETSC_FALSE,&m,
4727 (
const PetscInt**)&jj,&done); CHKERRQ(ierr);
4728 MFEM_VERIFY(done,
"Unable to perform MatRestoreRowIJ on "
4729 << i <<
" l2l matrix");
4732 for (
int i = 0; i < l2l.Size(); i++) { nl += l2l[i]->Width(); }
4734 const int* vdata = mark->
GetData();
4735 int* sdata = sub_dof_marker.
GetData();
4736 int cumh = 0, cumw = 0;
4737 for (
int i = 0; i < l2l.Size(); i++)
4742 l2l[i]->BooleanMultTranspose(vf_marker,sf_marker);
4743 cumh += l2l[i]->Height();
4744 cumw += l2l[i]->Width();
4746 ierr = Convert_Array_IS(comm,
false,&sub_dof_marker,st,is); CCHKERRQ(comm,ierr);
4747 for (
int i = 0; i < pl2l.
Size(); i++)
4751 PetscFunctionReturn(0);
4754 #if !defined(PETSC_HAVE_HYPRE)
4756 #if defined(HYPRE_MIXEDINT)
4757 #error "HYPRE_MIXEDINT not supported"
4760 #include "_hypre_parcsr_mv.h"
4761 static PetscErrorCode MatConvert_hypreParCSR_AIJ(hypre_ParCSRMatrix* hA,
Mat* pA)
4763 MPI_Comm comm = hypre_ParCSRMatrixComm(hA);
4764 hypre_CSRMatrix *hdiag,*hoffd;
4766 PetscInt *dii,*djj,*oii,*ojj,*iptr;
4769 PetscErrorCode ierr;
4771 PetscFunctionBeginUser;
4772 hdiag = hypre_ParCSRMatrixDiag(hA);
4773 hoffd = hypre_ParCSRMatrixOffd(hA);
4774 m = hypre_CSRMatrixNumRows(hdiag);
4775 n = hypre_CSRMatrixNumCols(hdiag);
4776 dnnz = hypre_CSRMatrixNumNonzeros(hdiag);
4777 onnz = hypre_CSRMatrixNumNonzeros(hoffd);
4778 ierr = PetscMalloc1(m+1,&dii); CHKERRQ(ierr);
4779 ierr = PetscMalloc1(dnnz,&djj); CHKERRQ(ierr);
4780 ierr = PetscMalloc1(dnnz,&da); CHKERRQ(ierr);
4781 ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*
sizeof(
PetscInt));
4783 ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*
sizeof(
PetscInt));
4785 ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*
sizeof(
PetscScalar));
4792 ierr = PetscSortIntWithScalarArray(nc,iptr,aptr); CHKERRQ(ierr);
4796 ierr = MPI_Comm_size(comm,&size); CHKERRQ(ierr);
4801 ierr = PetscMalloc1(m+1,&oii); CHKERRQ(ierr);
4802 ierr = PetscMalloc1(onnz,&ojj); CHKERRQ(ierr);
4803 ierr = PetscMalloc1(onnz,&oa); CHKERRQ(ierr);
4804 ierr = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*
sizeof(
PetscInt));
4806 offdj = hypre_CSRMatrixJ(hoffd);
4807 coffd = hypre_ParCSRMatrixColMapOffd(hA);
4808 for (i=0; i<onnz; i++) { ojj[i] = coffd[offdj[i]]; }
4809 ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*
sizeof(
PetscScalar));
4816 ierr = PetscSortIntWithScalarArray(nc,iptr,aptr); CHKERRQ(ierr);
4820 ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,
4821 djj,da,oii,ojj,oa,pA); CHKERRQ(ierr);
4827 ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,pA); CHKERRQ(ierr);
4833 void *ptrs[6] = {dii,djj,da,oii,ojj,oa};
4834 const char *names[6] = {
"_mfem_csr_dii",
4845 ierr = PetscContainerCreate(comm,&c); CHKERRQ(ierr);
4846 ierr = PetscContainerSetPointer(c,ptrs[i]); CHKERRQ(ierr);
4847 ierr = PetscContainerSetUserDestroy(c,__mfem_array_container_destroy);
4851 ierr = PetscContainerDestroy(&c); CHKERRQ(ierr);
4853 PetscFunctionReturn(0);
4856 static PetscErrorCode MatConvert_hypreParCSR_IS(hypre_ParCSRMatrix* hA,
Mat* pA)
4859 ISLocalToGlobalMapping rl2g,cl2g;
4861 hypre_CSRMatrix *hdiag,*hoffd;
4862 MPI_Comm comm = hypre_ParCSRMatrixComm(hA);
4864 const char *names[2] = {
"_mfem_csr_aux",
4868 PetscInt *col_map_offd,*hdi,*hdj,*hoi,*hoj;
4870 PetscInt cum,dr,dc,oc,str,stc,nnz,i,jd,jo;
4871 PetscErrorCode ierr;
4873 PetscFunctionBeginUser;
4875 str = hypre_ParCSRMatrixFirstRowIndex(hA);
4876 stc = hypre_ParCSRMatrixFirstColDiag(hA);
4877 hdiag = hypre_ParCSRMatrixDiag(hA);
4878 hoffd = hypre_ParCSRMatrixOffd(hA);
4879 dr = hypre_CSRMatrixNumRows(hdiag);
4880 dc = hypre_CSRMatrixNumCols(hdiag);
4881 nnz = hypre_CSRMatrixNumNonzeros(hdiag);
4882 hdi = hypre_CSRMatrixI(hdiag);
4883 hdj = hypre_CSRMatrixJ(hdiag);
4884 hdd = hypre_CSRMatrixData(hdiag);
4885 oc = hypre_CSRMatrixNumCols(hoffd);
4886 nnz += hypre_CSRMatrixNumNonzeros(hoffd);
4887 hoi = hypre_CSRMatrixI(hoffd);
4888 hoj = hypre_CSRMatrixJ(hoffd);
4889 hod = hypre_CSRMatrixData(hoffd);
4892 ierr = ISCreateStride(comm,dr,str,1,&is); CHKERRQ(ierr);
4893 ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g); CHKERRQ(ierr);
4894 ierr = ISDestroy(&is); CHKERRQ(ierr);
4895 col_map_offd = hypre_ParCSRMatrixColMapOffd(hA);
4896 ierr = PetscMalloc1(dc+oc,&aux); CHKERRQ(ierr);
4897 for (i=0; i<dc; i++) { aux[i] = i+stc; }
4898 for (i=0; i<oc; i++) { aux[i+dc] = col_map_offd[i]; }
4899 ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is); CHKERRQ(ierr);
4900 ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g); CHKERRQ(ierr);
4901 ierr = ISDestroy(&is); CHKERRQ(ierr);
4904 ierr = MatCreate(comm,pA); CHKERRQ(ierr);
4905 ierr = MatSetSizes(*pA,dr,dc,PETSC_DECIDE,PETSC_DECIDE); CHKERRQ(ierr);
4906 ierr = MatSetType(*pA,MATIS); CHKERRQ(ierr);
4907 ierr = MatSetLocalToGlobalMapping(*pA,rl2g,cl2g); CHKERRQ(ierr);
4908 ierr = ISLocalToGlobalMappingDestroy(&rl2g); CHKERRQ(ierr);
4909 ierr = ISLocalToGlobalMappingDestroy(&cl2g); CHKERRQ(ierr);
4912 ierr = PetscMalloc1(nnz+dr+1,&aux); CHKERRQ(ierr);
4913 ierr = PetscMalloc1(nnz,&data); CHKERRQ(ierr);
4917 *ii = *(hdi++) + *(hoi++);
4918 for (jd=0,jo=0,cum=0; *ii<nnz; cum++)
4922 for (; jd<*hdi; jd++) { *jj++ = *hdj++; *aa++ = *hdd++; }
4923 for (; jo<*hoi; jo++) { *jj++ = *hoj++ + dc; *aa++ = *hod++; }
4924 *(++ii) = *(hdi++) + *(hoi++);
4925 ierr = PetscSortIntWithScalarArray(jd+jo-nc,jold,aold); CHKERRQ(ierr);
4927 for (; cum<dr; cum++) { *(++ii) = nnz; }
4931 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,ii,jj,aa,&lA);
4939 ierr = PetscContainerCreate(PETSC_COMM_SELF,&c); CHKERRQ(ierr);
4940 ierr = PetscContainerSetPointer(c,ptrs[i]); CHKERRQ(ierr);
4941 ierr = PetscContainerSetUserDestroy(c,__mfem_array_container_destroy);
4945 ierr = PetscContainerDestroy(&c); CHKERRQ(ierr);
4947 ierr = MatISSetLocalMat(*pA,lA); CHKERRQ(ierr);
4948 ierr = MatDestroy(&lA); CHKERRQ(ierr);
4949 ierr = MatAssemblyBegin(*pA,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
4950 ierr = MatAssemblyEnd(*pA,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
4951 PetscFunctionReturn(0);
4955 #endif // MFEM_USE_PETSC
4956 #endif // MFEM_USE_MPI
int Size() const
Logical size of the array.
void SetPreconditioner(Solver &precond)
void SetPrintLevel(int plev)
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
PetscParVector * B
Right-hand side and solution vector.
void CreatePrivateContext()
void Print(const char *fname=NULL, bool binary=false) const
Prints the matrix (to stdout if fname is NULL)
Vector()
Default constructor for Vector. Sets size = 0 and data = NULL.
MPI_Comm GetComm() const
MPI communicator.
bool isImplicit() const
True if type is IMPLICIT or HOMOGENEOUS.
bool isHomogeneous() const
True if type is HOMOGENEOUS.
virtual void SetOperator(const Operator &op)
Specification of the nonlinear operator.
void MakeDataOwner() const
Set the Vector data (host pointer) ownership flag.
Abstract class for PETSc's preconditioners.
PetscParMatrix & operator+=(const PetscParMatrix &B)
virtual void SetOperator(const Operator &op)
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
bool clcustom
Boolean to handle SetFromOptions calls.
void SetObjective(void(*obj)(Operator *op, const Vector &x, double *f))
Specification of an objective function to be used for line search.
void ConvertOperator(MPI_Comm comm, const Operator &op, Mat *B, Operator::Type tid)
virtual void MultTranspose(const Vector &b, Vector &x) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Wrapper for PETSc's matrix class.
PetscParVector & operator+=(const PetscParVector &y)
void ApplyBC(const Vector &x, Vector &y)
y = x on ess_tdof_list_c and y = g (internally evaluated) on ess_tdof_list
Abstract class for PETSc's linear solvers.
Vec x
The actual PETSc object.
void SetSize(int s)
Resize the vector to size s.
Abstract class for PETSc's solvers.
PetscODESolver::Type GetType() const
virtual void Run(Vector &x, double &t, double &dt, double t_final)
Perform time integration from time t [in] to time tf [in].
void SetPostCheck(void(*post)(Operator *op, const Vector &X, Vector &Y, Vector &W, bool &changed_y, bool &changed_w))
Base abstract class for first order time dependent operators.
PetscInt GetNumCols() const
Returns the local number of columns.
virtual Operator & GetGradient(const Vector &x) const
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate ...
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
void SetType(PetscODESolver::Type)
virtual void MonitorSolver(PetscSolver *solver)
Generic monitor to take access to the solver.
Pointer to an Operator of a specified type.
int * GetJ()
Return the array J.
void Mult(double a, const Vector &x, double b, Vector &y) const
Matvec: y = a A x + b y.
PetscBDDCSolver(MPI_Comm comm, Operator &op, const PetscBDDCSolverParams &opts=PetscBDDCSolverParams(), const std::string &prefix=std::string())
PetscParVector & AddValues(const Array< PetscInt > &, const Array< PetscScalar > &)
Add values in a vector.
void SetMat(Mat newA)
Replace the inner Mat Object. The reference count of newA is increased.
virtual ~PetscSolver()
Destroy the PetscParVectors allocated (if any).
PetscNonlinearSolver(MPI_Comm comm, const std::string &prefix=std::string())
T * GetData()
Returns the data.
virtual void SetTime(const double _t)
Set the current time.
virtual Operator & GetExplicitGradient(const Vector &x) const
Return an Operator representing dG/dx at the given point x and the currently set time.
int Size() const
Returns the size of the vector.
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
int * GetI()
Return the array I.
void Init()
Initialize with defaults. Does not initialize inherited members.
PetscPreconditioner(MPI_Comm comm, const std::string &prefix=std::string())
PetscClassId cid
The class id of the actual PETSc object.
ID for class PetscParMatrix, MATHYPRE format.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
int IsZeroBlock(int i, int j) const
Check if block (i,j) is a zero block.
Abstract parallel finite element space.
void SetRelTol(double tol)
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
void Randomize(PetscInt seed=0)
Set random values.
PetscPCGSolver(MPI_Comm comm, const std::string &prefix=std::string())
void Destroy()
Delete all owned data. Does not perform re-initialization with defaults.
void SetMonitor(PetscSolverMonitor *ctx)
Sets user-defined monitoring routine.
void SetMaxIter(int max_iter)
double * GetData() const
Return a pointer to the beginning of the Vector data.
PetscParVector & operator*=(PetscScalar d)
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
void SortColumnIndices()
Sort the column indices corresponding to each row.
void EliminateRowsCols(const Array< int > &rows_cols, const PetscParVector &X, PetscParVector &B, double diag=1.)
Eliminate rows and columns from the matrix, and rows from the vector B. Modify B with the BC values i...
Wrapper for PETSc's vector class.
void Customize(bool customize=true) const
Customize object with options set.
PetscInt M() const
Returns the global number of rows.
Operator & GetBlock(int i, int j)
Return a reference to block i,j.
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
virtual void Mult(const Vector &b, Vector &x) const
Application of the solver.
virtual ~PetscLinearSolver()
double * GetData()
Return the element data, i.e. the array A.
PetscInt GetNumRows() const
Returns the local number of rows.
virtual void ExplicitMult(const Vector &x, Vector &y) const
Perform the action of the explicit part of the operator, G: y = G(x, t) where t is the current time...
void Print(const char *fname=NULL, bool binary=false) const
Prints the vector (to stdout if fname is NULL)
void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B)
void ResetArray()
Reset the PETSc Vec object to use its default data. Call this method after the use of PlaceArray()...
PetscParVector & SetValues(const Array< PetscInt > &, const Array< PetscScalar > &)
Set values in a vector.
virtual void Init(TimeDependentOperator &f_, enum PetscODESolver::Type type)
Initialize the ODE solver.
PetscParMatrix * Transpose(bool action=false)
Returns the transpose of the PetscParMatrix.
ID for the base class Operator, i.e. any type.
PetscObject obj
The actual PETSc object (KSP, PC, SNES or TS).
Abstract class for monitoring PETSc's solvers.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
void SetJacobianType(Operator::Type type)
int to_int(const std::string &str)
Auxiliary class for BDDC customization.
virtual ~PetscPreconditioner()
ID for class PetscParMatrix, unspecified format.
MPI_Comm GetComm() const
Get the associated MPI communicator.
void FixResidualBC(const Vector &x, Vector &y)
y = x-g on ess_tdof_list, the rest of y is unchanged
void LoseData()
NULL-ifies the data.
PetscParMatrix & operator=(const PetscParMatrix &B)
PetscInt GetColStart() const
Returns the global index of the first local column.
virtual ~PetscParVector()
Calls PETSc's destroy function.
PetscFieldSplitSolver(MPI_Comm comm, Operator &op, const std::string &prefix=std::string())
PetscBCHandler * bchandler
Handler for boundary conditions.
PetscParMatrix & operator-=(const PetscParMatrix &B)
virtual Solver * NewPreconditioner(const OperatorHandle &oh)=0
void ScaleRows(const Vector &s)
Scale the local row i by s(i).
void Zero(Vector &x)
Replace boundary dofs with 0.
void Assign(const T *)
Copy data from a pointer. Size() elements are copied.
MPI_Comm GetComm() const
Get the associated MPI communicator.
void Shift(double s)
Shift diagonal by a constant.
ID for class PetscParMatrix, MATNEST format.
void SetTime(double t)
Sets the current time.
ID for class PetscParMatrix, MATSHELL format.
void trans(const Vector &x, Vector &p)
PetscParMatrix()
Create an empty matrix to be used as a reference to an existing matrix.
void SetBCHandler(PetscBCHandler *bch)
Sets the object to handle essential boundary conditions.
Wrapper for hypre's parallel vector class.
void ZeroBC(const Vector &x, Vector &y)
y = x on ess_tdof_list_c and y = 0 on ess_tdof_list
virtual void Eval(double t, Vector &g)
Boundary conditions evaluation.
Type
Enumeration defining IDs for some classes derived from Operator.
PetscBCHandler(Type _type=ZERO)
PetscInt NNZ() const
Returns the number of nonzeros.
PetscLinearSolver(MPI_Comm comm, const std::string &prefix=std::string(), bool wrap=true)
virtual void Mult(const Vector &b, Vector &x) const
Application of the preconditioner.
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Mat A
The actual PETSc object.
virtual Operator & GetImplicitGradient(const Vector &x, const Vector &k, double shift) const
Return an Operator representing (dF/dk shift + dF/dx) at the given x, k, and the currently set time...
PetscInt N() const
Returns the global number of columns.
PetscParVector * GetY() const
Returns the inner vector in the range of A (it creates it if needed)
void SetJacobianType(Operator::Type type)
void MakeWrapper(MPI_Comm comm, const Operator *op, Mat *B)
Creates a wrapper around a mfem::Operator op using PETSc's MATSHELL object and returns the Mat in B...
void SetUpdate(void(*update)(Operator *op, int it, const mfem::Vector &F, const mfem::Vector &X, const mfem::Vector &D, const mfem::Vector &P))
bool operatorset
Boolean to handle SetOperator calls.
PetscParVector * X
Auxiliary vectors for typecasting.
void MultTranspose(double a, const Vector &x, double b, Vector &y) const
Matvec transpose: y = a A^T x + b y.
int height
Dimension of the output / number of rows in the matrix.
PetscParVector & operator-=(const PetscParVector &y)
virtual void ImplicitMult(const Vector &x, const Vector &k, Vector &y) const
Perform the action of the implicit part of the operator, F: y = F(x, k, t) where t is the current tim...
virtual ~PetscNonlinearSolver()
virtual ~PetscODESolver()
PetscSolver()
Construct an empty PetscSolver. Initialize protected objects to NULL.
void SetPreconditionerFactory(PetscPreconditionerFactory *factory)
Sets the object for the creation of the preconditioner.
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
Mat ReleaseMat(bool dereference)
Release the PETSc Mat object. If dereference is true, decrement the refcount of the Mat object...
PetscParVector(MPI_Comm comm, PetscInt glob_size, PetscInt *col=NULL)
Creates vector with given global size and partitioning of the columns.
virtual void MonitorSolution(PetscInt it, PetscReal norm, const Vector &x)
Monitor the solution vector x.
ID for class PetscParMatrix, MATAIJ format.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
MPI_Comm GetComm() const
Get the associated MPI communicator.
PetscInt GlobalSize() const
Returns the global number of rows.
PetscInt GetRowStart() const
Returns the global index of the first local row.
struct _p_PetscObject * PetscObject
void SetAbsTol(double tol)
Helper class for handling essential boundary conditions.
PetscParVector * GetX() const
Returns the inner vector in the domain of A (it creates it if needed)
Array< int > & GetTDofs()
Gets essential dofs (local, per-process numbering)
PetscParVector & operator=(PetscScalar d)
Set constant values.
HYPRE_Int * GetTrueDofOffsets() const
virtual void MonitorResidual(PetscInt it, PetscReal norm, const Vector &r)
Monitor the residual vector r.
int NumRowBlocks() const
Return the number of row blocks.
PetscParMatrix * TripleMatrixProduct(PetscParMatrix *R, PetscParMatrix *A, PetscParMatrix *P)
Returns the matrix R * A * P.
Identity Operator I: x -> x.
void PlaceArray(PetscScalar *temp_data)
Temporarily replace the data of the PETSc Vec object. To return to the original data array...
void FreePrivateContext()
TimeDependentOperator * f
Pointer to the associated TimeDependentOperator.
Vector * GlobalVector() const
Returns the global vector in each processor.
PetscODESolver(MPI_Comm comm, const std::string &prefix=std::string())
void EliminateRows(const Array< int > &rows)
Eliminate only the rows from the matrix.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Wrapper for hypre's ParCSR matrix class.
A class to handle Block systems in a matrix-free implementation.
virtual void Mult(const Vector &b, Vector &x) const
Application of the solver.
void operator*=(double s)
Scale all entries by s: A_scaled = s*A.
int NumColBlocks() const
Return the number of column blocks.
void SetTDofs(Array< int > &list)
Sets essential dofs (local, per-process numbering)
ID for class PetscParMatrix, MATIS format.
virtual void MultTranspose(const Vector &b, Vector &x) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
int width
Dimension of the input / number of columns in the matrix.
void * private_ctx
Private context for solver.
void SetUp(PetscInt n)
SetUp the helper object, where n is the size of the solution vector.
void MakeRef(const PetscParMatrix &master)
Makes this object a reference to another PetscParMatrix.
void ScaleCols(const Vector &s)
Scale the local col i by s(i).