14 #include "../config/config.hpp"
20 #include "../fem/fem.hpp"
21 #if defined(PETSC_HAVE_HYPRE)
22 #include "petscmathypre.h"
34 #define PCHKERRQ(obj,err) do { \
37 PetscError(PetscObjectComm((PetscObject)(obj)),__LINE__,_MFEM_FUNC_NAME, \
38 __FILE__,(err),PETSC_ERROR_REPEAT,NULL); \
39 MFEM_ABORT("Error in PETSc. See stacktrace above."); \
42 #define CCHKERRQ(comm,err) do { \
45 PetscError(comm,__LINE__,_MFEM_FUNC_NAME, \
46 __FILE__,(err),PETSC_ERROR_REPEAT,NULL); \
47 MFEM_ABORT("Error in PETSc. See stacktrace above."); \
52 static PetscErrorCode __mfem_ts_monitor(TS,PetscInt,PetscReal,Vec,
void*);
53 static PetscErrorCode __mfem_ts_rhsfunction(TS,PetscReal,Vec,Vec,
void*);
54 static PetscErrorCode __mfem_ts_rhsjacobian(TS,PetscReal,Vec,Mat,Mat,
56 static PetscErrorCode __mfem_ts_ifunction(TS,PetscReal,Vec,Vec,Vec,
void*);
57 static PetscErrorCode __mfem_ts_ijacobian(TS,PetscReal,Vec,Vec,
60 static PetscErrorCode __mfem_snes_monitor(SNES,PetscInt,PetscReal,
void*);
61 static PetscErrorCode __mfem_snes_jacobian(SNES,Vec,Mat,Mat,
void*);
62 static PetscErrorCode __mfem_snes_function(SNES,Vec,Vec,
void*);
63 static PetscErrorCode __mfem_ksp_monitor(KSP,PetscInt,PetscReal,
void*);
64 static PetscErrorCode __mfem_pc_shell_apply(PC,Vec,Vec);
65 static PetscErrorCode __mfem_pc_shell_apply_transpose(PC,Vec,Vec);
66 static PetscErrorCode __mfem_pc_shell_setup(PC);
67 static PetscErrorCode __mfem_pc_shell_destroy(PC);
68 static PetscErrorCode __mfem_pc_shell_view(PC,PetscViewer);
69 static PetscErrorCode __mfem_mat_shell_apply(Mat,Vec,Vec);
70 static PetscErrorCode __mfem_mat_shell_apply_transpose(Mat,Vec,Vec);
71 static PetscErrorCode __mfem_mat_shell_destroy(Mat);
72 static PetscErrorCode __mfem_array_container_destroy(
void*);
73 static PetscErrorCode __mfem_matarray_container_destroy(
void*);
76 static PetscErrorCode Convert_Array_IS(MPI_Comm,
bool,
const mfem::Array<int>*,
80 static PetscErrorCode MakeShellPC(PC,
mfem::Solver&,
bool);
81 static PetscErrorCode MakeShellPCWithFactory(PC,
87 #if !defined(PETSC_HAVE_HYPRE)
88 static PetscErrorCode MatConvert_hypreParCSR_AIJ(hypre_ParCSRMatrix*,Mat*);
89 static PetscErrorCode MatConvert_hypreParCSR_IS(hypre_ParCSRMatrix*,Mat*);
96 } __mfem_mat_shell_ctx;
103 unsigned long int numprec;
104 } __mfem_pc_shell_ctx;
121 PetscReal cached_shift;
122 PetscObjectState cached_ijacstate;
123 PetscObjectState cached_rhsjacstate;
127 static PetscErrorCode ierr;
136 void PetscParVector::_SetDataAndSize_()
138 const PetscScalar *array;
141 ierr = VecGetArrayRead(x,&array); PCHKERRQ(x,ierr);
142 ierr = VecGetLocalSize(x,&n); PCHKERRQ(x,ierr);
143 SetDataAndSize((PetscScalar*)array,n);
144 ierr = VecRestoreArrayRead(x,&array); PCHKERRQ(x,ierr);
147 PetscInt PetscParVector::GlobalSize()
const
150 ierr = VecGetSize(x,&N); PCHKERRQ(x,ierr);
154 PetscParVector::PetscParVector(MPI_Comm comm,
const Vector &_x) :
Vector()
156 ierr = VecCreate(comm,&
x); CCHKERRQ(comm,ierr);
157 ierr = VecSetSizes(
x,_x.
Size(),PETSC_DECIDE); PCHKERRQ(
x,ierr);
158 ierr = VecSetType(
x,VECSTANDARD); PCHKERRQ(
x,ierr);
165 ierr = VecCreate(comm,&
x); CCHKERRQ(comm,ierr);
169 MPI_Comm_rank(comm, &myid);
170 ierr = VecSetSizes(
x,col[myid+1]-col[myid],PETSC_DECIDE); PCHKERRQ(
x,ierr);
174 ierr = VecSetSizes(
x,PETSC_DECIDE,glob_size); PCHKERRQ(
x,ierr);
176 ierr = VecSetType(
x,VECSTANDARD); PCHKERRQ(
x,ierr);
182 MPI_Comm comm = PetscObjectComm((PetscObject)
x);
183 ierr = VecDestroy(&x); CCHKERRQ(comm,ierr);
187 PetscScalar *_data, PetscInt *col) :
Vector()
189 MFEM_VERIFY(col,
"Missing distribution");
191 MPI_Comm_rank(comm, &myid);
192 ierr = VecCreateMPIWithArray(comm,1,col[myid+1]-col[myid],glob_size,_data,
193 &
x); CCHKERRQ(comm,ierr)
199 ierr = VecDuplicate(y.
x,&
x); PCHKERRQ(
x,ierr);
204 bool transpose,
bool allocate) :
Vector()
206 PetscInt loc = transpose ? op.
Height() : op.
Width();
209 ierr = VecCreate(comm,&
x);
211 ierr = VecSetSizes(
x,loc,PETSC_DECIDE);
213 ierr = VecSetType(
x,VECSTANDARD);
220 ierr = VecCreateMPIWithArray(comm,1,loc,PETSC_DECIDE,NULL,
221 &
x); CCHKERRQ(comm,ierr);
227 bool transpose,
bool allocate) :
Vector()
232 ierr = MatCreateVecs(pA,&
x,NULL); PCHKERRQ(pA,ierr);
236 ierr = MatCreateVecs(pA,NULL,&
x); PCHKERRQ(pA,ierr);
240 ierr = VecReplaceArray(
x,NULL); PCHKERRQ(
x,ierr);
249 ierr = PetscObjectReference((PetscObject)y); PCHKERRQ(y,ierr);
259 MPI_Comm comm = pfes->
GetComm();
260 ierr = VecCreate(comm,&
x); CCHKERRQ(comm,ierr);
262 PetscMPIInt myid = 0;
263 if (!HYPRE_AssumedPartitionCheck())
265 MPI_Comm_rank(comm,&myid);
267 ierr = VecSetSizes(
x,offsets[myid+1]-offsets[myid],PETSC_DECIDE);
269 ierr = VecSetType(
x,VECSTANDARD); PCHKERRQ(
x,ierr);
280 ierr = VecScatterCreateToAll(
x,&scctx,&vout); PCHKERRQ(
x,ierr);
281 ierr = VecScatterBegin(scctx,
x,vout,INSERT_VALUES,SCATTER_FORWARD);
283 ierr = VecScatterEnd(scctx,
x,vout,INSERT_VALUES,SCATTER_FORWARD);
285 ierr = VecScatterDestroy(&scctx); PCHKERRQ(
x,ierr);
286 ierr = VecGetArray(vout,&array); PCHKERRQ(
x,ierr);
287 ierr = VecGetLocalSize(vout,&size); PCHKERRQ(
x,ierr);
290 ierr = VecRestoreArray(vout,&array); PCHKERRQ(
x,ierr);
291 ierr = VecDestroy(&vout); PCHKERRQ(
x,ierr);
300 ierr = VecSet(
x,d); PCHKERRQ(
x,ierr);
306 ierr = VecCopy(y.
x,
x); PCHKERRQ(
x,ierr);
312 ierr = VecPlaceArray(
x,temp_data); PCHKERRQ(
x,ierr);
317 ierr = VecResetArray(
x); PCHKERRQ(
x,ierr);
324 ierr = PetscRandomCreate(PetscObjectComm((PetscObject)
x),&rctx);
326 ierr = PetscRandomSetSeed(rctx,(
unsigned long)seed); PCHKERRQ(x,ierr);
327 ierr = PetscRandomSeed(rctx); PCHKERRQ(x,ierr);
328 ierr = VecSetRandom(x,rctx); PCHKERRQ(x,ierr);
329 ierr = PetscRandomDestroy(&rctx); PCHKERRQ(x,ierr);
340 ierr = PetscViewerBinaryOpen(PetscObjectComm((PetscObject)
x),fname,
341 FILE_MODE_WRITE,&view);
345 ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)
x),fname,&view);
348 ierr = VecView(
x,view); PCHKERRQ(
x,ierr);
349 ierr = PetscViewerDestroy(&view); PCHKERRQ(
x,ierr);
353 ierr = VecView(
x,NULL); PCHKERRQ(
x,ierr);
362 ierr = MatGetLocalSize(
A,&N,NULL); PCHKERRQ(
A,ierr);
369 ierr = MatGetLocalSize(
A,NULL,&N); PCHKERRQ(
A,ierr);
376 ierr = MatGetSize(
A,&N,NULL); PCHKERRQ(
A,ierr);
383 ierr = MatGetSize(
A,NULL,&N); PCHKERRQ(
A,ierr);
390 ierr = MatGetInfo(
A,MAT_GLOBAL_SUM,&info); PCHKERRQ(
A,ierr);
391 return (PetscInt)info.nz_used;
428 BlockDiagonalConstructor(comm,row_starts,row_starts,diag,
436 PetscInt global_num_cols, PetscInt *row_starts,
441 BlockDiagonalConstructor(comm,row_starts,col_starts,diag,
452 MPI_Comm comm = PetscObjectComm((PetscObject)
A);
453 ierr = MatDestroy(&A); CCHKERRQ(comm,ierr);
460 #if defined(PETSC_HAVE_HYPRE)
461 ierr = MatCreateFromParCSR(B,MATAIJ,PETSC_USE_POINTER,&
A);
463 ierr = MatConvert_hypreParCSR_AIJ(B,&
A); CCHKERRQ(B.
GetComm(),ierr);
472 MPI_Comm comm = PetscObjectComm((PetscObject)
A);
473 ierr = MatDestroy(&A); CCHKERRQ(comm,ierr);
480 ierr = MatDuplicate(B,MAT_COPY_VALUES,&
A); CCHKERRQ(B.
GetComm(),ierr);
488 ierr = MatDuplicate(B,MAT_COPY_VALUES,&
A); CCHKERRQ(B.
GetComm(),ierr);
492 MFEM_VERIFY(
height == B.
Height(),
"Invalid number of local rows");
493 MFEM_VERIFY(
width == B.
Width(),
"Invalid number of local columns");
494 ierr = MatAXPY(
A,1.0,B,DIFFERENT_NONZERO_PATTERN); CCHKERRQ(B.
GetComm(),ierr);
499 void PetscParMatrix::
500 BlockDiagonalConstructor(MPI_Comm comm,
501 PetscInt *row_starts, PetscInt *col_starts,
505 PetscInt lrsize,lcsize,rstart,cstart;
506 PetscMPIInt myid = 0,commsize;
508 ierr = MPI_Comm_size(comm,&commsize); CCHKERRQ(comm,ierr);
509 if (!HYPRE_AssumedPartitionCheck())
511 ierr = MPI_Comm_rank(comm,&myid); CCHKERRQ(comm,ierr);
513 lrsize = row_starts[myid+1]-row_starts[myid];
514 rstart = row_starts[myid];
515 lcsize = col_starts[myid+1]-col_starts[myid];
516 cstart = col_starts[myid];
521 ierr = ISCreateStride(comm,diag->
Height(),rstart,1,&is); CCHKERRQ(comm,ierr);
522 ISLocalToGlobalMapping rl2g,cl2g;
523 ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g); PCHKERRQ(is,ierr);
524 ierr = ISDestroy(&is); CCHKERRQ(comm,ierr);
525 if (row_starts != col_starts)
527 ierr = ISCreateStride(comm,diag->
Width(),cstart,1,&is);
529 ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g); PCHKERRQ(is,ierr);
530 ierr = ISDestroy(&is); CCHKERRQ(comm,ierr);
534 ierr = PetscObjectReference((PetscObject)rl2g); PCHKERRQ(rl2g,ierr);
539 ierr = MatCreate(comm,&A); CCHKERRQ(comm,ierr);
540 ierr = MatSetSizes(A,lrsize,lcsize,PETSC_DECIDE,PETSC_DECIDE);
542 ierr = MatSetType(A,MATIS); PCHKERRQ(A,ierr);
543 ierr = MatSetLocalToGlobalMapping(A,rl2g,cl2g); PCHKERRQ(A,ierr);
544 ierr = ISLocalToGlobalMappingDestroy(&rl2g); PCHKERRQ(A,ierr)
545 ierr = ISLocalToGlobalMappingDestroy(&cl2g); PCHKERRQ(A,ierr)
549 ierr = MatISGetLocalMat(A,&lA); PCHKERRQ(A,ierr);
550 if (sizeof(PetscInt) == sizeof(
int))
552 ierr = MatSeqAIJSetPreallocationCSR(lA,diag->
GetI(),diag->
GetJ(),
553 diag->
GetData()); PCHKERRQ(lA,ierr);
557 MFEM_ABORT(
"64bit indices not yet supported");
563 PetscInt *dii,*djj,*oii,
569 ierr = PetscMalloc1(m,&dii); CCHKERRQ(PETSC_COMM_SELF,ierr);
570 ierr = PetscMalloc1(nnz,&djj); CCHKERRQ(PETSC_COMM_SELF,ierr);
571 ierr = PetscMalloc1(nnz,&da); CCHKERRQ(PETSC_COMM_SELF,ierr);
572 if (
sizeof(PetscInt) ==
sizeof(
int))
574 ierr = PetscMemcpy(dii,diag->
GetI(),m*
sizeof(PetscInt));
575 CCHKERRQ(PETSC_COMM_SELF,ierr);
576 ierr = PetscMemcpy(djj,diag->
GetJ(),nnz*
sizeof(PetscInt));
577 CCHKERRQ(PETSC_COMM_SELF,ierr);
578 ierr = PetscMemcpy(da,diag->
GetData(),nnz*
sizeof(PetscScalar));
579 CCHKERRQ(PETSC_COMM_SELF,ierr);
583 MFEM_ABORT(
"64bit indices not yet supported");
585 ierr = PetscCalloc1(m,&oii);
586 CCHKERRQ(PETSC_COMM_SELF,ierr);
589 ierr = MatCreateMPIAIJWithSplitArrays(comm,lrsize,lcsize,PETSC_DECIDE,
591 dii,djj,da,oii,NULL,NULL,&A);
596 ierr = MatCreateSeqAIJWithArrays(comm,lrsize,lcsize,dii,djj,da,&A);
600 void *ptrs[4] = {dii,djj,da,oii};
601 const char *names[4] = {
"_mfem_csr_dii",
606 for (PetscInt i=0; i<4; i++)
610 ierr = PetscContainerCreate(comm,&c); CCHKERRQ(comm,ierr);
611 ierr = PetscContainerSetPointer(c,ptrs[i]); CCHKERRQ(comm,ierr);
612 ierr = PetscContainerSetUserDestroy(c,__mfem_array_container_destroy);
614 ierr = PetscObjectCompose((PetscObject)A,names[i],(PetscObject)c);
616 ierr = PetscContainerDestroy(&c); CCHKERRQ(comm,ierr);
621 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); PCHKERRQ(A,ierr);
622 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); PCHKERRQ(A,ierr);
637 __mfem_mat_shell_ctx *ctx =
new __mfem_mat_shell_ctx;
638 ierr = MatCreate(comm,A); CCHKERRQ(comm,ierr);
640 PETSC_DECIDE,PETSC_DECIDE); PCHKERRQ(A,ierr);
641 ierr = MatSetType(*A,MATSHELL); PCHKERRQ(A,ierr);
642 ierr = MatShellSetContext(*A,(
void *)ctx); PCHKERRQ(A,ierr);
643 ierr = MatShellSetOperation(*A,MATOP_MULT,
644 (
void (*)())__mfem_mat_shell_apply);
646 ierr = MatShellSetOperation(*A,MATOP_MULT_TRANSPOSE,
647 (
void (*)())__mfem_mat_shell_apply_transpose);
649 ierr = MatShellSetOperation(*A,MATOP_DESTROY,
650 (
void (*)())__mfem_mat_shell_destroy);
652 ierr = MatSetUp(*A); PCHKERRQ(*A,ierr);
653 ctx->op =
const_cast<Operator *
>(op);
671 PetscBool ismatis,istrans;
673 ierr = PetscObjectTypeCompare((PetscObject)(pA->
A),MATTRANSPOSEMAT,&istrans);
679 ierr = PetscObjectReference((PetscObject)(pA->
A));
684 ierr = PetscObjectTypeCompare((PetscObject)(pA->
A),MATIS,&ismatis);
689 ierr = MatTransposeGetMat(pA->
A,&At); CCHKERRQ(pA->
GetComm(),ierr);
690 ierr = PetscObjectTypeCompare((PetscObject)(At),MATIS,&ismatis);
703 ierr = MatISGetMPIXAIJ(At,MAT_INITIAL_MATRIX,&B); PCHKERRQ(pA->
A,ierr);
704 ierr = MatCreateTranspose(B,A); PCHKERRQ(pA->
A,ierr);
705 ierr = MatDestroy(&B); PCHKERRQ(pA->
A,ierr);
709 ierr = MatISGetMPIXAIJ(pA->
A,MAT_INITIAL_MATRIX,A);
710 PCHKERRQ(pA->
A,ierr);
716 ierr = MPI_Comm_size(comm,&size); CCHKERRQ(comm,ierr);
722 ierr = MatConvert(At,size > 1 ? MATMPIAIJ : MATSEQAIJ,MAT_INITIAL_MATRIX,&B);
723 PCHKERRQ(pA->
A,ierr);
724 ierr = MatCreateTranspose(B,A); PCHKERRQ(pA->
A,ierr);
725 ierr = MatDestroy(&B); PCHKERRQ(pA->
A,ierr);
729 ierr = MatConvert(pA->
A, size > 1 ? MATMPIAIJ : MATSEQAIJ,MAT_INITIAL_MATRIX,A);
730 PCHKERRQ(pA->
A,ierr);
739 ierr = MatConvert(At,MATIS,MAT_INITIAL_MATRIX,&B); PCHKERRQ(pA->
A,ierr);
740 ierr = MatCreateTranspose(B,A); PCHKERRQ(pA->
A,ierr);
741 ierr = MatDestroy(&B); PCHKERRQ(pA->
A,ierr);
745 ierr = MatConvert(pA->
A,MATIS,MAT_INITIAL_MATRIX,A); PCHKERRQ(pA->
A,ierr);
748 #if defined(PETSC_HAVE_HYPRE)
754 ierr = MatConvert(At,MATHYPRE,MAT_INITIAL_MATRIX,&B); PCHKERRQ(pA->
A,ierr);
755 ierr = MatCreateTranspose(B,A); PCHKERRQ(pA->
A,ierr);
756 ierr = MatDestroy(&B); PCHKERRQ(pA->
A,ierr);
760 ierr = MatConvert(pA->
A,MATHYPRE,MAT_INITIAL_MATRIX,A); PCHKERRQ(pA->
A,ierr);
770 MFEM_ABORT(
"Unsupported operator type conversion " << tid)
777 #if defined(PETSC_HAVE_HYPRE)
778 ierr = MatCreateFromParCSR(const_cast<HypreParMatrix&>(*pH),MATAIJ,
779 PETSC_USE_POINTER,A);
781 ierr = MatConvert_hypreParCSR_AIJ(const_cast<HypreParMatrix&>(*pH),A);
787 #if defined(PETSC_HAVE_HYPRE)
788 ierr = MatCreateFromParCSR(const_cast<HypreParMatrix&>(*pH),MATIS,
789 PETSC_USE_POINTER,A);
791 ierr = MatConvert_hypreParCSR_IS(const_cast<HypreParMatrix&>(*pH),A);
795 #if defined(PETSC_HAVE_HYPRE)
798 ierr = MatCreateFromParCSR(const_cast<HypreParMatrix&>(*pH),MATHYPRE,
799 PETSC_USE_POINTER,A);
809 MFEM_ABORT(
"Conversion from HypreParCSR to operator type = " << tid <<
810 " is not implemented");
815 Mat *mats,*matsl2l = NULL;
820 ierr = PetscCalloc1(nr*nc,&mats); CCHKERRQ(PETSC_COMM_SELF,ierr);
823 ierr = PetscCalloc1(nr,&matsl2l); CCHKERRQ(PETSC_COMM_SELF,ierr);
827 PetscBool needl2l = PETSC_TRUE;
837 ierr = PetscObjectQuery((PetscObject)mats[i*nc+j],
"__mfem_l2l",
839 PCHKERRQ(mats[i*nc+j],ierr);
847 ierr = PetscContainerGetPointer(c,(
void**)&l2l);
849 MFEM_VERIFY(l2l->
Size() == 1,
"Unexpected size "
850 << l2l->
Size() <<
" for block row " << i );
851 ierr = PetscObjectReference((PetscObject)(*l2l)[0]);
853 matsl2l[i] = (*l2l)[0];
854 needl2l = PETSC_FALSE;
860 ierr = MatCreateNest(comm,nr,NULL,nc,NULL,mats,A); CCHKERRQ(comm,ierr);
863 ierr = MatConvert(*A,MATIS,MAT_INPLACE_MATRIX,A); CCHKERRQ(comm,ierr);
866 for (PetscInt i=0; i<nr; i++) { (*vmatsl2l)[i] = matsl2l[i]; }
867 ierr = PetscFree(matsl2l); CCHKERRQ(PETSC_COMM_SELF,ierr);
870 ierr = PetscContainerCreate(comm,&c); CCHKERRQ(comm,ierr);
871 ierr = PetscContainerSetPointer(c,vmatsl2l); PCHKERRQ(c,ierr);
872 ierr = PetscContainerSetUserDestroy(c,__mfem_matarray_container_destroy);
874 ierr = PetscObjectCompose((PetscObject)(*A),
"__mfem_l2l",(PetscObject)c);
876 ierr = PetscContainerDestroy(&c); CCHKERRQ(comm,ierr);
878 for (i=0; i<nr*nc; i++) { ierr = MatDestroy(&mats[i]); CCHKERRQ(comm,ierr); }
879 ierr = PetscFree(mats); CCHKERRQ(PETSC_COMM_SELF,ierr);
883 MFEM_VERIFY(tid ==
PETSC_MATAIJ,
"Unsupported operation");
886 ierr = MatCreate(comm,A); CCHKERRQ(comm,ierr);
887 ierr = MatSetSizes(*A,pI->
Height(),pI->
Width(),PETSC_DECIDE,PETSC_DECIDE);
889 ierr = MatSetType(*A,MATAIJ); PCHKERRQ(*A,ierr);
890 ierr = MatMPIAIJSetPreallocation(*A,1,NULL,0,NULL); PCHKERRQ(*A,ierr);
891 ierr = MatSeqAIJSetPreallocation(*A,1,NULL); PCHKERRQ(*A,ierr);
892 ierr = MatSetOption(*A,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE); PCHKERRQ(*A,ierr);
893 ierr = MatGetOwnershipRange(*A,&rst,NULL); PCHKERRQ(*A,ierr);
894 for (PetscInt i = rst; i < rst+pI->
Height(); i++)
896 ierr = MatSetValue(*A,i,i,1.,INSERT_VALUES); PCHKERRQ(*A,ierr);
898 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY); PCHKERRQ(*A,ierr);
899 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY); PCHKERRQ(*A,ierr);
911 MPI_Comm comm = MPI_COMM_NULL;
912 ierr = PetscObjectGetComm((PetscObject)A,&comm); PCHKERRQ(A,ierr);
913 ierr = MatDestroy(&A); CCHKERRQ(comm,ierr);
924 ierr = PetscObjectReference((PetscObject)a); PCHKERRQ(a,ierr);
934 static void MatMultKernel(Mat A,PetscScalar a,Vec X,PetscScalar b,Vec Y,
937 PetscErrorCode (*f)(Mat,Vec,Vec);
938 PetscErrorCode (*fadd)(Mat,Vec,Vec,Vec);
941 f = MatMultTranspose;
942 fadd = MatMultTransposeAdd;
953 ierr = VecScale(X,a); PCHKERRQ(A,ierr);
954 ierr = (*fadd)(A,X,Y,Y); PCHKERRQ(A,ierr);
955 ierr = VecScale(X,1./a); PCHKERRQ(A,ierr);
959 ierr = VecScale(X,a); PCHKERRQ(A,ierr);
960 ierr = VecScale(Y,b); PCHKERRQ(A,ierr);
961 ierr = (*fadd)(A,X,Y,Y); PCHKERRQ(A,ierr);
962 ierr = VecScale(X,1./a); PCHKERRQ(A,ierr);
966 ierr = (*f)(A,X,Y); PCHKERRQ(A,ierr);
969 ierr = VecScale(Y,a); PCHKERRQ(A,ierr);
981 ierr = VecScale(Y,b); PCHKERRQ(A,ierr);
985 ierr = VecSet(Y,0.); PCHKERRQ(A,ierr);
992 ierr = PetscObjectReference((PetscObject)master.
A); PCHKERRQ(master.
A,ierr);
1004 MFEM_VERIFY(A,
"Mat not present");
1014 MFEM_VERIFY(A,
"Mat not present");
1025 ierr = MatCreateTranspose(A,&B); PCHKERRQ(A,ierr);
1029 ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&B); PCHKERRQ(A,ierr);
1036 ierr = MatScale(A,s); PCHKERRQ(A,ierr);
1041 MFEM_ASSERT(x.
Size() ==
Width(),
"invalid x.Size() = " << x.
Size()
1042 <<
", expected size = " <<
Width());
1043 MFEM_ASSERT(y.
Size() ==
Height(),
"invalid y.Size() = " << y.
Size()
1044 <<
", expected size = " <<
Height());
1050 MatMultKernel(A,a,XX->x,b,YY->
x,
false);
1058 MFEM_ASSERT(x.
Size() ==
Height(),
"invalid x.Size() = " << x.
Size()
1059 <<
", expected size = " <<
Height());
1060 MFEM_ASSERT(y.
Size() ==
Width(),
"invalid y.Size() = " << y.
Size()
1061 <<
", expected size = " <<
Width());
1067 MatMultKernel(A,a,YY->
x,b,XX->x,
true);
1080 ierr = PetscViewerBinaryOpen(PetscObjectComm((PetscObject)A),fname,
1081 FILE_MODE_WRITE,&view);
1085 ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)A),fname,&view);
1088 ierr = MatView(A,view); PCHKERRQ(A,ierr);
1089 ierr = PetscViewerDestroy(&view); PCHKERRQ(A,ierr);
1093 ierr = MatView(A,NULL); PCHKERRQ(A,ierr);
1100 Mat pA = *A,pP = *P,pRt = *Rt;
1102 PetscBool Aismatis,Pismatis,Rtismatis;
1105 "Petsc RAP: Number of local cols of A " << A->
Width() <<
1106 " differs from number of local rows of P " << P->
Height());
1108 "Petsc RAP: Number of local rows of A " << A->
Height() <<
1109 " differs from number of local rows of Rt " << Rt->
Height());
1110 ierr = PetscObjectTypeCompare((PetscObject)pA,MATIS,&Aismatis);
1112 ierr = PetscObjectTypeCompare((PetscObject)pP,MATIS,&Pismatis);
1114 ierr = PetscObjectTypeCompare((PetscObject)pRt,MATIS,&Rtismatis);
1121 ISLocalToGlobalMapping cl2gP,cl2gRt;
1122 PetscInt rlsize,clsize,rsize,csize;
1124 ierr = MatGetLocalToGlobalMapping(pP,NULL,&cl2gP); PCHKERRQ(pA,ierr);
1125 ierr = MatGetLocalToGlobalMapping(pRt,NULL,&cl2gRt); PCHKERRQ(pA,ierr);
1126 ierr = MatGetLocalSize(pP,NULL,&clsize); PCHKERRQ(pP,ierr);
1127 ierr = MatGetLocalSize(pRt,NULL,&rlsize); PCHKERRQ(pRt,ierr);
1128 ierr = MatGetSize(pP,NULL,&csize); PCHKERRQ(pP,ierr);
1129 ierr = MatGetSize(pRt,NULL,&rsize); PCHKERRQ(pRt,ierr);
1130 ierr = MatCreate(A->
GetComm(),&B); PCHKERRQ(pA,ierr);
1131 ierr = MatSetSizes(B,rlsize,clsize,rsize,csize); PCHKERRQ(B,ierr);
1132 ierr = MatSetType(B,MATIS); PCHKERRQ(B,ierr);
1133 ierr = MatSetLocalToGlobalMapping(B,cl2gRt,cl2gP); PCHKERRQ(B,ierr);
1134 ierr = MatISGetLocalMat(pA,&lA); PCHKERRQ(pA,ierr);
1135 ierr = MatISGetLocalMat(pP,&lP); PCHKERRQ(pA,ierr);
1136 ierr = MatISGetLocalMat(pRt,&lRt); PCHKERRQ(pA,ierr);
1139 ierr = MatPtAP(lA,lP,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&lB);
1145 ierr = MatTranspose(lRt,MAT_INITIAL_MATRIX,&lR); PCHKERRQ(lRt,ierr);
1146 ierr = MatMatMatMult(lR,lA,lP,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&lB);
1148 ierr = MatDestroy(&lR); PCHKERRQ(lRt,ierr);
1156 ierr = PetscObjectReference((PetscObject)lRt); PCHKERRQ(lRt,ierr);
1157 (*vmatsl2l)[0] = lRt;
1160 ierr = PetscContainerCreate(PetscObjectComm((PetscObject)B),&c);
1162 ierr = PetscContainerSetPointer(c,vmatsl2l); PCHKERRQ(c,ierr);
1163 ierr = PetscContainerSetUserDestroy(c,__mfem_matarray_container_destroy);
1165 ierr = PetscObjectCompose((PetscObject)B,
"__mfem_l2l",(PetscObject)c);
1167 ierr = PetscContainerDestroy(&c); PCHKERRQ(B,ierr);
1171 ierr = MatISSetLocalMat(B,lB); PCHKERRQ(lB,ierr);
1172 ierr = MatDestroy(&lB); PCHKERRQ(lA,ierr);
1173 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); PCHKERRQ(B,ierr);
1174 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); PCHKERRQ(B,ierr);
1180 ierr = MatPtAP(pA,pP,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B);
1186 ierr = MatTranspose(pRt,MAT_INITIAL_MATRIX,&pR); PCHKERRQ(Rt,ierr);
1187 ierr = MatMatMatMult(pR,pA,pP,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B);
1189 ierr = MatDestroy(&pR); PCHKERRQ(pRt,ierr);
1206 ierr = MatDuplicate(A,MAT_COPY_VALUES,&Ae); PCHKERRQ(A,ierr);
1208 ierr = MatAXPY(Ae,-1.,A,SAME_NONZERO_PATTERN); PCHKERRQ(A,ierr);
1216 MFEM_ABORT(
"PetscParMatrix::EliminateRowsCols(). To be implemented");
1224 ierr = MatGetSize(A,&M,&N); PCHKERRQ(A,ierr);
1225 MFEM_VERIFY(M == N,
"Rectangular case unsupported");
1228 ierr = MatSetOption(A,MAT_NO_OFF_PROC_ZERO_ROWS,PETSC_TRUE); PCHKERRQ(A,ierr);
1232 ierr = MatGetOwnershipRange(A,&rst,NULL); PCHKERRQ(A,ierr);
1235 ierr = Convert_Array_IS(
GetComm(),
true,&rows_cols,rst,&dir); PCHKERRQ(A,ierr);
1238 ierr = MatZeroRowsColumnsIS(A,dir,1.,NULL,NULL); PCHKERRQ(A,ierr);
1242 ierr = MatZeroRowsColumnsIS(A,dir,1.,X,B); PCHKERRQ(A,ierr);
1244 ierr = ISDestroy(&dir); PCHKERRQ(A,ierr);
1254 ierr = PetscObjectDereference((PetscObject)A); CCHKERRQ(comm,ierr);
1263 MFEM_VERIFY(A,
"no associated PETSc Mat object");
1264 PetscObject oA = (PetscObject)(this->A);
1266 ierr = PetscObjectTypeCompare(oA, MATAIJ, &ok); PCHKERRQ(A,ierr);
1268 ierr = PetscObjectTypeCompare(oA, MATSEQAIJ, &ok); PCHKERRQ(A,ierr);
1270 ierr = PetscObjectTypeCompare(oA, MATMPIAIJ, &ok); PCHKERRQ(A,ierr);
1272 ierr = PetscObjectTypeCompare(oA, MATIS, &ok); PCHKERRQ(A,ierr);
1274 ierr = PetscObjectTypeCompare(oA, MATSHELL, &ok); PCHKERRQ(A,ierr);
1276 ierr = PetscObjectTypeCompare(oA, MATNEST, &ok); PCHKERRQ(A,ierr);
1278 #if defined(PETSC_HAVE_HYPRE)
1279 ierr = PetscObjectTypeCompare(oA, MATHYPRE, &ok); PCHKERRQ(A,ierr);
1289 const PetscScalar *array;
1293 Ae.
Mult(-1.0, X, 1.0, B);
1296 ierr = MatGetDiagonal(pA,diag); PCHKERRQ(pA,ierr);
1297 ierr = VecGetArrayRead(diag,&array); PCHKERRQ(diag,ierr);
1298 for (
int i = 0; i < ess_dof_list.
Size(); i++)
1300 int r = ess_dof_list[i];
1301 B(r) = array[r] * X(r);
1303 ierr = VecRestoreArrayRead(diag,&array); PCHKERRQ(diag,ierr);
1332 if (
cid == KSP_CLASSID)
1335 ierr = KSPSetTolerances(ksp,tol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
1337 else if (
cid == SNES_CLASSID)
1339 SNES snes = (SNES)
obj;
1340 ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tol,PETSC_DEFAULT,PETSC_DEFAULT,
1343 else if (
cid == TS_CLASSID)
1346 ierr = TSSetTolerances(ts,PETSC_DECIDE,NULL,tol,NULL);
1350 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
1357 if (
cid == KSP_CLASSID)
1360 ierr = KSPSetTolerances(ksp,PETSC_DEFAULT,tol,PETSC_DEFAULT,PETSC_DEFAULT);
1362 else if (
cid == SNES_CLASSID)
1364 SNES snes = (SNES)
obj;
1365 ierr = SNESSetTolerances(snes,tol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,
1368 else if (
cid == TS_CLASSID)
1371 ierr = TSSetTolerances(ts,tol,NULL,PETSC_DECIDE,NULL);
1375 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
1382 if (
cid == KSP_CLASSID)
1385 ierr = KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,
1388 else if (
cid == SNES_CLASSID)
1390 SNES snes = (SNES)
obj;
1391 ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,
1392 max_iter,PETSC_DEFAULT);
1394 else if (
cid == TS_CLASSID)
1397 ierr = TSSetMaxSteps(ts,max_iter);
1401 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
1409 typedef PetscErrorCode (*myPetscFunc)(
void**);
1410 PetscViewerAndFormat *vf = NULL;
1411 PetscViewer viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm(
obj));
1415 ierr = PetscViewerAndFormatCreate(viewer,PETSC_VIEWER_DEFAULT,&vf);
1418 if (
cid == KSP_CLASSID)
1422 typedef PetscErrorCode (*myMonitor)(KSP,PetscInt,PetscReal,
void*);
1426 ierr = KSPMonitorCancel(ksp); PCHKERRQ(ksp,ierr);
1430 ierr = KSPMonitorSet(ksp,(myMonitor)KSPMonitorDefault,vf,
1431 (myPetscFunc)PetscViewerAndFormatDestroy);
1436 ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE); PCHKERRQ(ksp,ierr);
1437 ierr = KSPMonitorSet(ksp,(myMonitor)KSPMonitorSingularValue,vf,
1438 (myPetscFunc)PetscViewerAndFormatDestroy);
1442 ierr = PetscViewerAndFormatCreate(viewer,PETSC_VIEWER_DEFAULT,&vf);
1443 PCHKERRQ(viewer,ierr);
1444 ierr = KSPMonitorSet(ksp,(myMonitor)KSPMonitorTrueResidualNorm,vf,
1445 (myPetscFunc)PetscViewerAndFormatDestroy);
1450 else if (
cid == SNES_CLASSID)
1452 typedef PetscErrorCode (*myMonitor)(SNES,PetscInt,PetscReal,
void*);
1453 SNES snes = (SNES)
obj;
1456 ierr = SNESMonitorCancel(snes); PCHKERRQ(snes,ierr);
1460 ierr = SNESMonitorSet(snes,(myMonitor)SNESMonitorDefault,vf,
1461 (myPetscFunc)PetscViewerAndFormatDestroy);
1462 PCHKERRQ(snes,ierr);
1465 else if (
cid == TS_CLASSID)
1470 ierr = TSMonitorCancel(ts); PCHKERRQ(ts,ierr);
1475 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
1481 if (
cid == KSP_CLASSID)
1483 ierr = KSPMonitorSet((KSP)
obj,__mfem_ksp_monitor,ctx,NULL);
1486 else if (
cid == SNES_CLASSID)
1488 ierr = SNESMonitorSet((SNES)
obj,__mfem_snes_monitor,ctx,NULL);
1491 else if (
cid == TS_CLASSID)
1493 ierr = TSMonitorSet((TS)
obj,__mfem_ts_monitor,ctx,NULL);
1498 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
1505 if (
cid == SNES_CLASSID)
1507 __mfem_snes_ctx* snes_ctx = (__mfem_snes_ctx*)
private_ctx;
1510 else if (
cid == TS_CLASSID)
1512 __mfem_ts_ctx* ts_ctx = (__mfem_ts_ctx*)
private_ctx;
1517 MFEM_ABORT(
"Handling of essential bc only implemented for nonlinear and time-dependent solvers");
1524 if (
cid == TS_CLASSID)
1529 ierr = TSGetSNES((TS)
obj,&snes); PCHKERRQ(obj,ierr);
1530 ierr = SNESGetKSP(snes,&ksp); PCHKERRQ(obj,ierr);
1531 ierr = KSPGetPC(ksp,&pc); PCHKERRQ(obj,ierr);
1533 else if (
cid == SNES_CLASSID)
1537 ierr = SNESGetKSP((SNES)
obj,&ksp); PCHKERRQ(obj,ierr);
1538 ierr = KSPGetPC(ksp,&pc); PCHKERRQ(obj,ierr);
1540 else if (
cid == KSP_CLASSID)
1542 ierr = KSPGetPC((KSP)
obj,&pc); PCHKERRQ(obj,ierr);
1544 else if (
cid == PC_CLASSID)
1550 MFEM_ABORT(
"No support for PetscPreconditionerFactory for this object");
1552 ierr = MakeShellPCWithFactory(pc,factory); PCHKERRQ(pc,ierr);
1557 if (!customize) {
clcustom =
true; }
1560 if (
cid == PC_CLASSID)
1563 ierr = PCSetFromOptions(pc); PCHKERRQ(pc, ierr);
1565 else if (
cid == KSP_CLASSID)
1568 ierr = KSPSetFromOptions(ksp); PCHKERRQ(ksp, ierr);
1570 else if (
cid == SNES_CLASSID)
1572 SNES snes = (SNES)
obj;
1573 ierr = SNESSetFromOptions(snes); PCHKERRQ(snes, ierr);
1575 else if (
cid == TS_CLASSID)
1578 ierr = TSSetFromOptions(ts); PCHKERRQ(ts, ierr);
1582 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
1590 if (
cid == KSP_CLASSID)
1593 KSPConvergedReason reason;
1594 ierr = KSPGetConvergedReason(ksp,&reason);
1596 return reason > 0 ? 1 : 0;
1598 else if (
cid == SNES_CLASSID)
1600 SNES snes = (SNES)
obj;
1601 SNESConvergedReason reason;
1602 ierr = SNESGetConvergedReason(snes,&reason);
1603 PCHKERRQ(snes,ierr);
1604 return reason > 0 ? 1 : 0;
1606 else if (
cid == TS_CLASSID)
1609 TSConvergedReason reason;
1610 ierr = TSGetConvergedReason(ts,&reason);
1612 return reason > 0 ? 1 : 0;
1616 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
1623 if (
cid == KSP_CLASSID)
1627 ierr = KSPGetIterationNumber(ksp,&its);
1631 else if (
cid == SNES_CLASSID)
1633 SNES snes = (SNES)
obj;
1635 ierr = SNESGetIterationNumber(snes,&its);
1636 PCHKERRQ(snes,ierr);
1639 else if (
cid == TS_CLASSID)
1643 ierr = TSGetStepNumber(ts,&its);
1649 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
1656 if (
cid == KSP_CLASSID)
1660 ierr = KSPGetResidualNorm(ksp,&norm);
1664 if (
cid == SNES_CLASSID)
1666 SNES snes = (SNES)
obj;
1668 ierr = SNESGetFunctionNorm(snes,&norm);
1669 PCHKERRQ(snes,ierr);
1674 MFEM_ABORT(
"CLASSID = " <<
cid <<
" is not implemented!");
1675 return PETSC_MAX_REAL;
1682 if (
cid == SNES_CLASSID)
1684 __mfem_snes_ctx *snes_ctx;
1685 ierr = PetscNew(&snes_ctx); CCHKERRQ(PETSC_COMM_SELF,ierr);
1686 snes_ctx->op = NULL;
1687 snes_ctx->bchandler = NULL;
1688 snes_ctx->work = NULL;
1692 else if (
cid == TS_CLASSID)
1694 __mfem_ts_ctx *ts_ctx;
1695 ierr = PetscNew(&ts_ctx); CCHKERRQ(PETSC_COMM_SELF,ierr);
1697 ts_ctx->bchandler = NULL;
1698 ts_ctx->work = NULL;
1699 ts_ctx->cached_shift = std::numeric_limits<PetscReal>::min();
1700 ts_ctx->cached_ijacstate = -1;
1701 ts_ctx->cached_rhsjacstate = -1;
1712 if (
cid == SNES_CLASSID)
1714 __mfem_snes_ctx *snes_ctx = (__mfem_snes_ctx *)
private_ctx;
1715 delete snes_ctx->work;
1717 else if (
cid == TS_CLASSID)
1719 __mfem_ts_ctx *ts_ctx = (__mfem_ts_ctx *)
private_ctx;
1720 delete ts_ctx->work;
1722 ierr = PetscFree(
private_ctx); CCHKERRQ(PETSC_COMM_SELF,ierr);
1729 : bctype(_type), setup(false), eval_t(0.0),
1730 eval_t_cached(std::numeric_limits<double>::min())
1738 ess_tdof_list.
Assign(list);
1744 if (setup) {
return; }
1748 this->
Eval(eval_t,eval_g);
1749 eval_t_cached = eval_t;
1760 if (!setup) { MFEM_ABORT(
"PetscBCHandler not yet setup"); }
1764 for (PetscInt i = 0; i < ess_tdof_list.
Size(); ++i)
1766 y[ess_tdof_list[i]] = 0.0;
1771 if (bctype !=
CONSTANT && eval_t != eval_t_cached)
1773 Eval(eval_t,eval_g);
1774 eval_t_cached = eval_t;
1776 for (PetscInt i = 0; i < ess_tdof_list.
Size(); ++i)
1778 y[ess_tdof_list[i]] = eval_g[ess_tdof_list[i]];
1785 if (!setup) { MFEM_ABORT(
"PetscBCHandler not yet setup"); }
1788 for (PetscInt i = 0; i < ess_tdof_list.
Size(); ++i)
1790 y[ess_tdof_list[i]] = x[ess_tdof_list[i]];
1795 for (PetscInt i = 0; i < ess_tdof_list.
Size(); ++i)
1797 y[ess_tdof_list[i]] = x[ess_tdof_list[i]] - eval_g[ess_tdof_list[i]];
1809 ierr = KSPCreate(comm,&ksp); CCHKERRQ(comm,ierr);
1810 obj = (PetscObject)ksp;
1811 ierr = PetscObjectGetClassId(
obj,&
cid); PCHKERRQ(
obj,ierr);
1812 ierr = KSPSetOptionsPrefix(ksp, prefix.c_str()); PCHKERRQ(ksp, ierr);
1816 const std::string &prefix)
1821 obj = (PetscObject)ksp;
1822 ierr = PetscObjectGetClassId(obj,&
cid); PCHKERRQ(obj,ierr);
1823 ierr = KSPSetOptionsPrefix(ksp, prefix.c_str()); PCHKERRQ(ksp, ierr);
1828 const std::string &prefix)
1833 obj = (PetscObject)ksp;
1834 ierr = PetscObjectGetClassId(obj, &
cid); PCHKERRQ(obj, ierr);
1835 ierr = KSPSetOptionsPrefix(ksp, prefix.c_str()); PCHKERRQ(ksp, ierr);
1850 ierr = KSPGetOperatorsSet(ksp,NULL,&pmat); PCHKERRQ(ksp,ierr);
1853 ierr = KSPGetOperators(ksp,NULL,&P); PCHKERRQ(ksp,ierr);
1854 ierr = PetscObjectReference((PetscObject)P); PCHKERRQ(ksp,ierr);
1858 bool delete_pA =
false;
1876 MFEM_VERIFY(pA,
"Unsupported operation!");
1883 PetscInt nheight,nwidth,oheight,owidth;
1885 ierr = KSPGetOperators(ksp,&C,NULL); PCHKERRQ(ksp,ierr);
1886 ierr = MatGetSize(A,&nheight,&nwidth); PCHKERRQ(A,ierr);
1887 ierr = MatGetSize(C,&oheight,&owidth); PCHKERRQ(A,ierr);
1888 if (nheight != oheight || nwidth != owidth)
1892 ierr = KSPReset(ksp); PCHKERRQ(ksp,ierr);
1901 ierr = KSPSetOperators(ksp,A,P); PCHKERRQ(ksp,ierr);
1902 ierr = MatDestroy(&P); PCHKERRQ(ksp,ierr);
1906 ierr = KSPSetOperators(ksp,A,A); PCHKERRQ(ksp,ierr);
1916 if (delete_pA) {
delete pA; }
1931 bool delete_pA =
false;
1949 MFEM_VERIFY(pA,
"Unsupported operation!");
1952 bool delete_ppA =
false;
1955 if (oA == poA && !wrap)
1965 MFEM_VERIFY(ppA,
"Unsupported operation!");
1974 PetscInt nheight,nwidth,oheight,owidth;
1976 ierr = KSPGetOperators(ksp,&C,NULL); PCHKERRQ(ksp,ierr);
1977 ierr = MatGetSize(A,&nheight,&nwidth); PCHKERRQ(A,ierr);
1978 ierr = MatGetSize(C,&oheight,&owidth); PCHKERRQ(A,ierr);
1979 if (nheight != oheight || nwidth != owidth)
1983 ierr = KSPReset(ksp); PCHKERRQ(ksp,ierr);
1990 ierr = KSPSetOperators(ksp,A,P); PCHKERRQ(ksp,ierr);
1999 if (delete_pA) {
delete pA; }
2000 if (delete_ppA) {
delete ppA; }
2010 ierr = KSPGetOperatorsSet(ksp,&amat,NULL); PCHKERRQ(ksp,ierr);
2013 ierr = KSPGetOperators(ksp,&A,NULL); PCHKERRQ(ksp,ierr);
2014 ierr = PetscObjectReference((PetscObject)A); PCHKERRQ(ksp,ierr);
2019 ierr = KSPSetPC(ksp,*ppc); PCHKERRQ(ksp,ierr);
2028 ierr = KSPGetPC(ksp,&pc); PCHKERRQ(ksp,ierr);
2029 ierr = MakeShellPC(pc,precond,
false); PCHKERRQ(ksp,ierr);
2035 ierr = KSPGetOperators(ksp,NULL,&P); PCHKERRQ(ksp,ierr);
2036 ierr = PetscObjectReference((PetscObject)P); PCHKERRQ(ksp,ierr);
2037 ierr = KSPSetOperators(ksp,A,P); PCHKERRQ(ksp,ierr);
2038 ierr = MatDestroy(&A); PCHKERRQ(ksp,ierr);
2039 ierr = MatDestroy(&P); PCHKERRQ(ksp,ierr);
2050 ierr = KSPGetOperators(ksp, &pA, NULL); PCHKERRQ(
obj, ierr);
2067 ierr = KSPSetInitialGuessNonzero(ksp, (PetscBool)
iterative_mode);
2068 PCHKERRQ(ksp, ierr);
2071 ierr = KSPSolve(ksp,
B->
x, X->x); PCHKERRQ(ksp,ierr);
2080 ierr = PetscObjectGetComm((PetscObject)ksp,&comm); PCHKERRQ(ksp,ierr);
2081 ierr = KSPDestroy(&ksp); CCHKERRQ(comm,ierr);
2090 ierr = KSPSetType(ksp,KSPCG); PCHKERRQ(ksp,ierr);
2092 ierr = KSPSetNormType(ksp,KSP_NORM_NATURAL); PCHKERRQ(ksp,ierr);
2099 ierr = KSPSetType(ksp,KSPCG); PCHKERRQ(ksp,ierr);
2101 ierr = KSPSetNormType(ksp,KSP_NORM_NATURAL); PCHKERRQ(ksp,ierr);
2105 const std::string &prefix)
2109 ierr = KSPSetType(ksp,KSPCG); PCHKERRQ(ksp,ierr);
2111 ierr = KSPSetNormType(ksp,KSP_NORM_NATURAL); PCHKERRQ(ksp,ierr);
2117 const std::string &prefix)
2121 ierr = PCCreate(comm,&pc); CCHKERRQ(comm,ierr);
2122 obj = (PetscObject)pc;
2123 ierr = PetscObjectGetClassId(
obj,&
cid); PCHKERRQ(
obj,ierr);
2124 ierr = PCSetOptionsPrefix(pc, prefix.c_str()); PCHKERRQ(pc, ierr);
2128 const string &prefix)
2133 obj = (PetscObject)pc;
2134 ierr = PetscObjectGetClassId(obj,&
cid); PCHKERRQ(obj,ierr);
2135 ierr = PCSetOptionsPrefix(pc, prefix.c_str()); PCHKERRQ(pc, ierr);
2140 const string &prefix)
2144 ierr = PCCreate(comm,&pc); CCHKERRQ(comm,ierr);
2145 obj = (PetscObject)pc;
2146 ierr = PetscObjectGetClassId(
obj,&
cid); PCHKERRQ(
obj,ierr);
2147 ierr = PCSetOptionsPrefix(pc, prefix.c_str()); PCHKERRQ(pc, ierr);
2153 bool delete_pA =
false;
2170 PetscInt nheight,nwidth,oheight,owidth;
2172 ierr = PCGetOperators(pc,&C,NULL); PCHKERRQ(pc,ierr);
2173 ierr = MatGetSize(A,&nheight,&nwidth); PCHKERRQ(A,ierr);
2174 ierr = MatGetSize(C,&oheight,&owidth); PCHKERRQ(A,ierr);
2175 if (nheight != oheight || nwidth != owidth)
2179 ierr = PCReset(pc); PCHKERRQ(pc,ierr);
2185 ierr = PCSetOperators(pc,pA->
A,pA->
A); PCHKERRQ(
obj,ierr);
2194 if (delete_pA) {
delete pA; };
2204 ierr = PCGetOperators(pc, NULL, &pA); PCHKERRQ(
obj, ierr);
2222 ierr = PCApply(pc,
B->
x, X->x); PCHKERRQ(pc, ierr);
2231 ierr = PetscObjectGetComm((PetscObject)pc,&comm); PCHKERRQ(pc,ierr);
2232 ierr = PCDestroy(&pc); CCHKERRQ(comm,ierr);
2239 MPI_Comm comm = PetscObjectComm(
obj);
2244 ierr = PCGetOperators(pc,NULL,&pA); PCHKERRQ(pc,ierr);
2248 ierr = PetscObjectTypeCompare((PetscObject)pA,MATIS,&ismatis);
2250 MFEM_VERIFY(ismatis,
"PetscBDDCSolver needs the matrix in unassembled format");
2253 ierr = PCSetType(pc,PCBDDC); PCHKERRQ(
obj,ierr);
2260 IS dir = NULL, neu = NULL;
2269 ierr = PetscObjectQuery((PetscObject)pA,
"__mfem_l2l",(PetscObject*)&c);
2270 MFEM_VERIFY(c,
"Local-to-local PETSc container not present");
2271 ierr = PetscContainerGetPointer(c,(
void**)&l2l); PCHKERRQ(c,ierr);
2278 PetscBool lpr = PETSC_FALSE,pr;
2279 if (opts.
ess_dof) { lpr = PETSC_TRUE; }
2280 ierr = MPI_Allreduce(&lpr,&pr,1,MPIU_BOOL,MPI_LOR,comm);
2282 MFEM_VERIFY(lpr == pr,
"ess_dof should be collectively set");
2284 if (opts.
nat_dof) { lpr = PETSC_TRUE; }
2285 ierr = MPI_Allreduce(&lpr,&pr,1,MPIU_BOOL,MPI_LOR,comm);
2287 MFEM_VERIFY(lpr == pr,
"nat_dof should be collectively set");
2289 PetscInt ms[2],Ms[2];
2290 ms[0] = -nf; ms[1] = nf;
2291 ierr = MPI_Allreduce(&ms,&Ms,2,MPIU_INT,MPI_MAX,comm);
2293 MFEM_VERIFY(-Ms[0] == Ms[1],
2294 "number of fields should be the same across processes");
2299 ierr = MatGetOwnershipRange(pA,&rst,NULL); PCHKERRQ(pA,ierr);
2306 ierr = Convert_Array_IS(comm,
true,opts.
ess_dof,st,&dir);
2307 CCHKERRQ(comm,ierr);
2308 ierr = PCBDDCSetDirichletBoundaries(pc,dir); PCHKERRQ(pc,ierr);
2313 ierr = Convert_Vmarks_IS(comm,*l2l,opts.
ess_dof,st,&dir);
2314 CCHKERRQ(comm,ierr);
2315 ierr = PCBDDCSetDirichletBoundariesLocal(pc,dir); PCHKERRQ(pc,ierr);
2324 ierr = Convert_Array_IS(comm,
true,opts.
nat_dof,st,&neu);
2325 CCHKERRQ(comm,ierr);
2326 ierr = PCBDDCSetNeumannBoundaries(pc,neu); PCHKERRQ(pc,ierr);
2331 ierr = Convert_Vmarks_IS(comm,*l2l,opts.
nat_dof,st,&neu);
2332 CCHKERRQ(comm,ierr);
2333 ierr = PCBDDCSetNeumannBoundariesLocal(pc,neu); PCHKERRQ(pc,ierr);
2340 ierr = PCBDDCSetDofsSplitting(pc,nf,fields); PCHKERRQ(pc,ierr);
2341 for (
int i = 0; i < nf; i++)
2343 ierr = ISDestroy(&fields[i]); CCHKERRQ(comm,ierr);
2345 ierr = PetscFree(fields); PCHKERRQ(pc,ierr);
2353 ParFiniteElementSpace *fespace = opts.
fespace;
2356 const FiniteElementCollection *fec = fespace->
FEColl();
2357 bool edgespace, rtspace;
2358 bool needint =
false;
2359 bool tracespace, rt_tracespace, edge_tracespace;
2361 PetscBool B_is_Trans = PETSC_FALSE;
2363 ParMesh *pmesh = (ParMesh *) fespace->GetMesh();
2364 dim = pmesh->Dimension();
2367 rtspace =
dynamic_cast<const RT_FECollection*
>(fec);
2368 edgespace =
dynamic_cast<const ND_FECollection*
>(fec);
2369 edge_tracespace =
dynamic_cast<const ND_Trace_FECollection*
>(fec);
2370 rt_tracespace =
dynamic_cast<const RT_Trace_FECollection*
>(fec);
2371 tracespace = edge_tracespace || rt_tracespace;
2374 if (fespace->GetNE() > 0)
2378 p = fespace->GetOrder(0);
2382 p = fespace->GetFaceOrder(0);
2383 if (dim == 2) { p++; }
2394 MFEM_WARNING(
"Tracespace case doesn't work for H(curl) and p=2,"
2395 " not using auxiliary quadrature");
2401 FiniteElementCollection *vfec;
2404 vfec =
new H1_Trace_FECollection(p,dim);
2408 vfec =
new H1_FECollection(p,dim);
2410 ParFiniteElementSpace *vfespace =
new ParFiniteElementSpace(pmesh,vfec);
2411 ParDiscreteLinearOperator *grad;
2412 grad =
new ParDiscreteLinearOperator(vfespace,fespace);
2415 grad->AddTraceFaceInterpolator(
new GradientInterpolator);
2419 grad->AddDomainInterpolator(
new GradientInterpolator);
2423 HypreParMatrix *hG = grad->ParallelAssemble();
2424 PetscParMatrix *G =
new PetscParMatrix(hG,
PETSC_MATAIJ);
2428 PetscBool conforming = PETSC_TRUE;
2429 if (pmesh->Nonconforming()) { conforming = PETSC_FALSE; }
2430 ierr = PCBDDCSetDiscreteGradient(pc,*G,p,0,PETSC_TRUE,conforming);
2442 MFEM_WARNING(
"Tracespace case doesn't work for H(div), not using"
2443 " auxiliary quadrature");
2452 PetscParMatrix *
B = NULL;
2458 FiniteElementCollection *auxcoll;
2459 if (tracespace) { auxcoll =
new RT_Trace_FECollection(p,dim); }
2460 else { auxcoll =
new L2_FECollection(p,dim); };
2461 ParFiniteElementSpace *pspace =
new ParFiniteElementSpace(pmesh,auxcoll);
2462 ParMixedBilinearForm *b =
new ParMixedBilinearForm(fespace,pspace);
2468 b->AddTraceFaceIntegrator(
new VectorFECurlIntegrator);
2472 b->AddDomainIntegrator(
new VectorFECurlIntegrator);
2479 b->AddTraceFaceIntegrator(
new VectorFEDivergenceIntegrator);
2483 b->AddDomainIntegrator(
new VectorFEDivergenceIntegrator);
2489 b->ParallelAssemble(Bh);
2491 Bh.SetOperatorOwner(
false);
2496 ierr = MatTranspose(pB,MAT_INPLACE_MATRIX,&pB); PCHKERRQ(pA,ierr);
2499 ierr = MatZeroRowsIS(pB,dir,0.,NULL,NULL); PCHKERRQ(pA,ierr);
2503 ierr = MatZeroRowsLocalIS(pB,dir,0.,NULL,NULL); PCHKERRQ(pA,ierr);
2505 B_is_Trans = PETSC_TRUE;
2514 ierr = PCBDDCSetDivergenceMat(pc,*B,B_is_Trans,NULL); PCHKERRQ(pc,ierr);
2518 ierr = ISDestroy(&dir); PCHKERRQ(pc,ierr);
2519 ierr = ISDestroy(&neu); PCHKERRQ(pc,ierr);
2524 const std::string &prefix)
2527 BDDCSolverConstructor(opts);
2533 const std::string &prefix)
2536 BDDCSolverConstructor(opts);
2541 const string &prefix)
2547 ierr = PCGetOperators(pc,&pA,NULL); PCHKERRQ(pc,ierr);
2552 ierr = PetscObjectTypeCompare((PetscObject)pA,MATNEST,&isnest);
2555 "PetscFieldSplitSolver needs the matrix in nested format.");
2559 ierr = PCSetType(pc,PCFIELDSPLIT); PCHKERRQ(pc,ierr);
2560 ierr = MatNestGetSize(pA,&nr,NULL); PCHKERRQ(pc,ierr);
2561 ierr = PetscCalloc1(nr,&isrow); CCHKERRQ(PETSC_COMM_SELF,ierr);
2562 ierr = MatNestGetISs(pA,isrow,NULL); PCHKERRQ(pc,ierr);
2567 for (PetscInt i=0; i<nr; i++)
2569 ierr = PCFieldSplitSetIS(pc,NULL,isrow[i]); PCHKERRQ(pc,ierr);
2571 ierr = PetscFree(isrow); CCHKERRQ(PETSC_COMM_SELF,ierr);
2577 const std::string &prefix)
2582 ierr = SNESCreate(comm, &snes); CCHKERRQ(comm, ierr);
2583 obj = (PetscObject)snes;
2584 ierr = PetscObjectGetClassId(
obj, &
cid); PCHKERRQ(
obj, ierr);
2585 ierr = SNESSetOptionsPrefix(snes, prefix.c_str()); PCHKERRQ(snes, ierr);
2592 const std::string &prefix)
2597 ierr = SNESCreate(comm, &snes); CCHKERRQ(comm, ierr);
2598 obj = (PetscObject)snes;
2599 ierr = PetscObjectGetClassId(
obj, &
cid); PCHKERRQ(
obj, ierr);
2600 ierr = SNESSetOptionsPrefix(snes, prefix.c_str()); PCHKERRQ(snes, ierr);
2611 SNES snes = (SNES)
obj;
2612 ierr = PetscObjectGetComm(
obj,&comm); PCHKERRQ(
obj, ierr);
2613 ierr = SNESDestroy(&snes); CCHKERRQ(comm, ierr);
2618 SNES snes = (SNES)
obj;
2625 ierr = SNESGetFunction(snes, NULL, NULL, &fctx);
2626 PCHKERRQ(snes, ierr);
2627 ierr = SNESGetJacobian(snes, NULL, NULL, NULL, &jctx);
2628 PCHKERRQ(snes, ierr);
2631 (
void*)&op == fctx &&
2632 (
void*)&op == jctx);
2633 ierr = MPI_Allreduce(&ls,&gs,1,MPIU_BOOL,MPI_LAND,
2634 PetscObjectComm((PetscObject)snes));
2635 PCHKERRQ(snes,ierr);
2638 ierr = SNESReset(snes); PCHKERRQ(snes,ierr);
2645 __mfem_snes_ctx *snes_ctx = (__mfem_snes_ctx*)
private_ctx;
2647 ierr = SNESSetFunction(snes, NULL, __mfem_snes_function, (
void *)snes_ctx);
2648 PCHKERRQ(snes, ierr);
2649 ierr = SNESSetJacobian(snes, NULL, NULL, __mfem_snes_jacobian,
2651 PCHKERRQ(snes, ierr);
2663 __mfem_snes_ctx *snes_ctx = (__mfem_snes_ctx*)
private_ctx;
2664 snes_ctx->jacType = jacType;
2669 SNES snes = (SNES)
obj;
2671 bool b_nonempty = b.
Size();
2675 if (b_nonempty) { B->PlaceArray(b.
GetData()); }
2685 ierr = SNESSolve(snes, B->x, X->x); PCHKERRQ(snes, ierr);
2687 if (b_nonempty) { B->ResetArray(); }
2697 ierr = TSCreate(comm,&ts); CCHKERRQ(comm,ierr);
2698 obj = (PetscObject)ts;
2699 ierr = PetscObjectGetClassId(
obj,&
cid); PCHKERRQ(
obj,ierr);
2700 ierr = TSSetOptionsPrefix(ts, prefix.c_str()); PCHKERRQ(ts, ierr);
2706 ierr = TSSetMaxSteps(ts,PETSC_MAX_INT-1);
2708 ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
2711 ierr = TSGetAdapt(ts,&tsad);
2713 ierr = TSAdaptSetType(tsad,TSADAPTNONE);
2721 ierr = PetscObjectGetComm(
obj,&comm); PCHKERRQ(
obj,ierr);
2722 ierr = TSDestroy(&ts); CCHKERRQ(comm,ierr);
2730 __mfem_ts_ctx *ts_ctx = (__mfem_ts_ctx*)
private_ctx;
2734 void *fctx = NULL,*jctx = NULL,*rfctx = NULL,*rjctx = NULL;
2738 ierr = TSGetIFunction(ts, NULL, NULL, &fctx);
2740 ierr = TSGetIJacobian(ts, NULL, NULL, NULL, &jctx);
2745 ierr = TSGetRHSFunction(ts, NULL, NULL, &rfctx);
2747 ierr = TSGetRHSJacobian(ts, NULL, NULL, NULL, &rjctx);
2757 ls = (PetscBool)(ls && (
void*)&f_ == fctx && (
void*)&f_ == jctx);
2761 ls = (PetscBool)(ls && (
void*)&f_ == rfctx && (
void*)&f_ == rjctx);
2763 ierr = MPI_Allreduce(&ls,&gs,1,MPIU_BOOL,MPI_LAND,
2764 PetscObjectComm((PetscObject)ts));
2768 ierr = TSReset(ts); PCHKERRQ(ts,ierr);
2771 ts_ctx->cached_shift = std::numeric_limits<PetscReal>::min();
2772 ts_ctx->cached_ijacstate = -1;
2773 ts_ctx->cached_rhsjacstate = -1;
2782 ierr = TSSetIFunction(ts, NULL, __mfem_ts_ifunction, (
void *)ts_ctx);
2784 ierr = TSSetIJacobian(ts, NULL, NULL, __mfem_ts_ijacobian, (
void *)ts_ctx);
2786 ierr = TSSetEquationType(ts, TS_EQ_IMPLICIT);
2793 ierr = TSSetEquationType(ts, TS_EQ_EXPLICIT);
2796 ierr = TSSetRHSFunction(ts, NULL, __mfem_ts_rhsfunction, (
void *)ts_ctx);
2798 ierr = TSSetRHSJacobian(ts, NULL, NULL, __mfem_ts_rhsjacobian, (
void *)ts_ctx);
2803 ts_ctx->type = type;
2806 ierr = TSSetProblemType(ts, TS_LINEAR);
2811 ierr = TSSetProblemType(ts, TS_NONLINEAR);
2818 __mfem_ts_ctx *ts_ctx = (__mfem_ts_ctx*)
private_ctx;
2819 ts_ctx->jacType = jacType;
2826 ierr = TSSetTime(ts, t); PCHKERRQ(ts, ierr);
2827 ierr = TSSetTimeStep(ts, dt); PCHKERRQ(ts, ierr);
2837 ierr = TSSetSolution(ts, *X); PCHKERRQ(ts, ierr);
2838 ierr = TSStep(ts); PCHKERRQ(ts, ierr);
2844 ierr = TSGetTime(ts,&pt); PCHKERRQ(ts,ierr);
2845 dt = pt - (PetscReal)t;
2853 ierr = TSSetTime(ts, t); PCHKERRQ(ts, ierr);
2854 ierr = TSSetTimeStep(ts, dt); PCHKERRQ(ts, ierr);
2855 ierr = TSSetMaxTime(ts, t_final); PCHKERRQ(ts, ierr);
2856 ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);
2867 ierr = TSSolve(ts, X->x); PCHKERRQ(ts, ierr);
2872 ierr = TSGetTime(ts, &pt); PCHKERRQ(ts,ierr);
2874 ierr = TSGetTimeStep(ts,&pt); PCHKERRQ(ts,ierr);
2880 #include "petsc/private/petscimpl.h"
2883 static PetscErrorCode __mfem_ts_monitor(TS ts, PetscInt it, PetscReal t, Vec x,
2888 PetscFunctionBeginUser;
2891 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,
"No monitor context provided");
2900 SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,
2901 "Cannot monitor the residual with TS");
2903 PetscFunctionReturn(0);
2906 static PetscErrorCode __mfem_ksp_monitor(KSP ksp, PetscInt it, PetscReal res,
2911 PetscErrorCode ierr;
2913 PetscFunctionBeginUser;
2916 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,
"No monitor context provided");
2920 ierr = KSPBuildSolution(ksp,NULL,&x); CHKERRQ(ierr);
2926 ierr = KSPBuildResidual(ksp,NULL,NULL,&x); CHKERRQ(ierr);
2930 PetscFunctionReturn(0);
2933 static PetscErrorCode __mfem_ts_ifunction(TS ts, PetscReal t, Vec x, Vec xp,
2936 __mfem_ts_ctx* ts_ctx = (__mfem_ts_ctx*)ctx;
2937 PetscErrorCode ierr;
2939 PetscFunctionBeginUser;
2947 if (ts_ctx->bchandler)
2950 if (!ts_ctx->work) { ts_ctx->work =
new mfem::Vector(xx.Size()); }
2966 ierr = PetscObjectStateIncrease((PetscObject)f); CHKERRQ(ierr);
2967 PetscFunctionReturn(0);
2970 static PetscErrorCode __mfem_ts_rhsfunction(TS ts, PetscReal t, Vec x, Vec f,
2973 __mfem_ts_ctx* ts_ctx = (__mfem_ts_ctx*)ctx;
2974 PetscErrorCode ierr;
2976 PetscFunctionBeginUser;
2977 if (ts_ctx->bchandler) { MFEM_ABORT(
"RHS evaluation with bc not implemented"); }
2987 ierr = PetscObjectStateIncrease((PetscObject)f); CHKERRQ(ierr);
2988 PetscFunctionReturn(0);
2991 static PetscErrorCode __mfem_ts_ijacobian(TS ts, PetscReal t, Vec x,
2992 Vec xp, PetscReal shift, Mat A, Mat P,
2995 __mfem_ts_ctx* ts_ctx = (__mfem_ts_ctx*)ctx;
2998 PetscReal eps = 0.001;
3000 PetscObjectState state;
3001 PetscErrorCode ierr;
3003 PetscFunctionBeginUser;
3011 ierr = PetscObjectStateGet((PetscObject)A,&state); CHKERRQ(ierr);
3013 std::abs(ts_ctx->cached_shift/shift - 1.0) < eps &&
3014 state == ts_ctx->cached_ijacstate) { PetscFunctionReturn(0); }
3017 ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr);
3018 ierr = VecGetArrayRead(xp,(
const PetscScalar**)&array); CHKERRQ(ierr);
3020 ierr = VecRestoreArrayRead(xp,(
const PetscScalar**)&array); CHKERRQ(ierr);
3021 ierr = VecGetArrayRead(x,(
const PetscScalar**)&array); CHKERRQ(ierr);
3022 if (!ts_ctx->bchandler)
3029 if (!ts_ctx->work) { ts_ctx->work =
new mfem::Vector(n); }
3036 ierr = VecRestoreArrayRead(x,(
const PetscScalar**)&array); CHKERRQ(ierr);
3040 if (!ts_ctx->bchandler) {
delete xx; }
3041 ts_ctx->cached_shift = shift;
3044 bool delete_pA =
false;
3047 if (!pA || pA->
GetType() != ts_ctx->jacType)
3055 if (ts_ctx->bchandler)
3065 ierr = MatHeaderReplace(A,&B); CHKERRQ(ierr);
3066 if (delete_pA) {
delete pA; }
3069 ierr = PetscObjectStateGet((PetscObject)A,&ts_ctx->cached_ijacstate);
3071 PetscFunctionReturn(0);
3074 static PetscErrorCode __mfem_ts_rhsjacobian(TS ts, PetscReal t, Vec x,
3075 Mat A, Mat P,
void *ctx)
3077 __mfem_ts_ctx* ts_ctx = (__mfem_ts_ctx*)ctx;
3081 PetscObjectState state;
3082 PetscErrorCode ierr;
3084 PetscFunctionBeginUser;
3090 ierr = PetscObjectStateGet((PetscObject)A,&state); CHKERRQ(ierr);
3092 state == ts_ctx->cached_rhsjacstate) { PetscFunctionReturn(0); }
3095 ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr);
3096 ierr = VecGetArrayRead(x,(
const PetscScalar**)&array); CHKERRQ(ierr);
3097 if (!ts_ctx->bchandler)
3104 if (!ts_ctx->work) { ts_ctx->work =
new mfem::Vector(n); }
3111 ierr = VecRestoreArrayRead(x,(
const PetscScalar**)&array); CHKERRQ(ierr);
3115 if (!ts_ctx->bchandler) {
delete xx; }
3118 bool delete_pA =
false;
3121 if (!pA || pA->
GetType() != ts_ctx->jacType)
3129 if (ts_ctx->bchandler)
3139 ierr = MatHeaderReplace(A,&B); CHKERRQ(ierr);
3140 if (delete_pA) {
delete pA; }
3145 ierr = TSRHSJacobianSetReuse(ts,PETSC_TRUE); PCHKERRQ(ts,ierr);
3147 ierr = PetscObjectStateGet((PetscObject)A,&ts_ctx->cached_rhsjacstate);
3149 PetscFunctionReturn(0);
3152 static PetscErrorCode __mfem_snes_monitor(SNES snes, PetscInt it, PetscReal res,
3157 PetscErrorCode ierr;
3159 PetscFunctionBeginUser;
3162 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,
"No monitor context provided");
3166 ierr = SNESGetSolution(snes,&x); CHKERRQ(ierr);
3172 ierr = SNESGetFunction(snes,&x,NULL,NULL); CHKERRQ(ierr);
3176 PetscFunctionReturn(0);
3179 static PetscErrorCode __mfem_snes_jacobian(SNES snes, Vec x, Mat A, Mat P,
3184 PetscErrorCode ierr;
3186 __mfem_snes_ctx *snes_ctx = (__mfem_snes_ctx*)ctx;
3188 PetscFunctionBeginUser;
3189 ierr = VecGetArrayRead(x,(
const PetscScalar**)&array); CHKERRQ(ierr);
3190 ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr);
3191 if (!snes_ctx->bchandler)
3198 if (!snes_ctx->work) { snes_ctx->work =
new mfem::Vector(n); }
3201 xx = snes_ctx->work;
3207 ierr = VecRestoreArrayRead(x,(
const PetscScalar**)&array); CHKERRQ(ierr);
3208 if (!snes_ctx->bchandler) {
delete xx; }
3211 bool delete_pA =
false;
3214 if (!pA || pA->
GetType() != snes_ctx->jacType)
3222 if (snes_ctx->bchandler)
3231 ierr = MatHeaderReplace(A,&B); CHKERRQ(ierr);
3232 if (delete_pA) {
delete pA; }
3233 PetscFunctionReturn(0);
3236 static PetscErrorCode __mfem_snes_function(SNES snes, Vec x, Vec f,
void *ctx)
3238 __mfem_snes_ctx* snes_ctx = (__mfem_snes_ctx*)ctx;
3240 PetscFunctionBeginUser;
3243 if (snes_ctx->bchandler)
3250 snes_ctx->op->Mult(*txx,ff);
3257 snes_ctx->op->Mult(xx,ff);
3260 ierr = PetscObjectStateIncrease((PetscObject)f); CHKERRQ(ierr);
3261 PetscFunctionReturn(0);
3264 static PetscErrorCode __mfem_mat_shell_apply(Mat A, Vec x, Vec y)
3266 __mfem_mat_shell_ctx *ctx;
3267 PetscErrorCode ierr;
3269 PetscFunctionBeginUser;
3270 ierr = MatShellGetContext(A,(
void **)&ctx); CHKERRQ(ierr);
3273 ctx->op->Mult(xx,yy);
3275 ierr = PetscObjectStateIncrease((PetscObject)y); CHKERRQ(ierr);
3276 PetscFunctionReturn(0);
3279 static PetscErrorCode __mfem_mat_shell_apply_transpose(Mat A, Vec x, Vec y)
3281 __mfem_mat_shell_ctx *ctx;
3282 PetscErrorCode ierr;
3284 PetscFunctionBeginUser;
3285 ierr = MatShellGetContext(A,(
void **)&ctx); CHKERRQ(ierr);
3288 ctx->op->MultTranspose(xx,yy);
3290 ierr = PetscObjectStateIncrease((PetscObject)y); CHKERRQ(ierr);
3291 PetscFunctionReturn(0);
3294 static PetscErrorCode __mfem_mat_shell_destroy(Mat A)
3296 __mfem_mat_shell_ctx *ctx;
3297 PetscErrorCode ierr;
3299 PetscFunctionBeginUser;
3300 ierr = MatShellGetContext(A,(
void **)&ctx); CHKERRQ(ierr);
3302 PetscFunctionReturn(0);
3305 static PetscErrorCode __mfem_pc_shell_view(PC pc, PetscViewer viewer)
3307 __mfem_pc_shell_ctx *ctx;
3308 PetscErrorCode ierr;
3310 PetscFunctionBeginUser;
3311 ierr = PCShellGetContext(pc,(
void **)&ctx); CHKERRQ(ierr);
3315 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
3322 ierr = PCView(*ppc,viewer); CHKERRQ(ierr);
3328 ierr = PetscViewerASCIIPrintf(viewer,
3329 "No information available on the mfem::Solver\n");
3333 if (isascii && ctx->factory)
3335 ierr = PetscViewerASCIIPrintf(viewer,
3336 "Number of preconditioners created by the factory %lu\n",ctx->numprec);
3340 PetscFunctionReturn(0);
3343 static PetscErrorCode __mfem_pc_shell_apply(PC pc, Vec x, Vec y)
3345 __mfem_pc_shell_ctx *ctx;
3346 PetscErrorCode ierr;
3348 PetscFunctionBeginUser;
3351 ierr = PCShellGetContext(pc,(
void **)&ctx); CHKERRQ(ierr);
3354 ctx->op->Mult(xx,yy);
3356 ierr = PetscObjectStateIncrease((PetscObject)y); CHKERRQ(ierr);
3362 PetscFunctionReturn(0);
3365 static PetscErrorCode __mfem_pc_shell_apply_transpose(PC pc, Vec x, Vec y)
3367 __mfem_pc_shell_ctx *ctx;
3368 PetscErrorCode ierr;
3370 PetscFunctionBeginUser;
3373 ierr = PCShellGetContext(pc,(
void **)&ctx); CHKERRQ(ierr);
3376 ctx->op->MultTranspose(xx,yy);
3378 ierr = PetscObjectStateIncrease((PetscObject)y); CHKERRQ(ierr);
3384 PetscFunctionReturn(0);
3387 static PetscErrorCode __mfem_pc_shell_setup(PC pc)
3389 __mfem_pc_shell_ctx *ctx;
3391 PetscFunctionBeginUser;
3392 ierr = PCShellGetContext(pc,(
void **)&ctx); CHKERRQ(ierr);
3403 ierr = PCGetOperators(pc,NULL,&B); CHKERRQ(ierr);
3412 PetscFunctionReturn(0);
3415 static PetscErrorCode __mfem_pc_shell_destroy(PC pc)
3417 __mfem_pc_shell_ctx *ctx;
3418 PetscErrorCode ierr;
3420 PetscFunctionBeginUser;
3421 ierr = PCShellGetContext(pc,(
void **)&ctx); CHKERRQ(ierr);
3427 PetscFunctionReturn(0);
3430 static PetscErrorCode __mfem_array_container_destroy(
void *ptr)
3432 PetscErrorCode ierr;
3434 PetscFunctionBeginUser;
3435 ierr = PetscFree(ptr); CHKERRQ(ierr);
3436 PetscFunctionReturn(0);
3439 static PetscErrorCode __mfem_matarray_container_destroy(
void *ptr)
3442 PetscErrorCode ierr;
3444 PetscFunctionBeginUser;
3445 for (
int i=0; i<a->
Size(); i++)
3448 MPI_Comm comm = PetscObjectComm((PetscObject)M);
3449 ierr = MatDestroy(&M); CCHKERRQ(comm,ierr);
3452 PetscFunctionReturn(0);
3457 PetscErrorCode MakeShellPC(PC pc,
mfem::Solver &precond,
bool ownsop)
3459 PetscFunctionBeginUser;
3460 __mfem_pc_shell_ctx *ctx =
new __mfem_pc_shell_ctx;
3462 ctx->ownsop = ownsop;
3463 ctx->factory = NULL;
3466 ierr = PCSetType(pc,PCSHELL); CHKERRQ(ierr);
3467 ierr = PCShellSetName(pc,
"MFEM Solver"); CHKERRQ(ierr);
3468 ierr = PCShellSetContext(pc,(
void *)ctx); CHKERRQ(ierr);
3469 ierr = PCShellSetApply(pc,__mfem_pc_shell_apply); CHKERRQ(ierr);
3470 ierr = PCShellSetApplyTranspose(pc,__mfem_pc_shell_apply_transpose);
3472 ierr = PCShellSetSetUp(pc,__mfem_pc_shell_setup); CHKERRQ(ierr);
3473 ierr = PCShellSetView(pc,__mfem_pc_shell_view); CHKERRQ(ierr);
3474 ierr = PCShellSetDestroy(pc,__mfem_pc_shell_destroy); CHKERRQ(ierr);
3475 PetscFunctionReturn(0);
3480 PetscErrorCode MakeShellPCWithFactory(PC pc,
3483 PetscFunctionBeginUser;
3484 __mfem_pc_shell_ctx *ctx =
new __mfem_pc_shell_ctx;
3487 ctx->factory = factory;
3490 ierr = PCSetType(pc,PCSHELL); CHKERRQ(ierr);
3491 ierr = PCShellSetName(pc,factory->
GetName()); CHKERRQ(ierr);
3492 ierr = PCShellSetContext(pc,(
void *)ctx); CHKERRQ(ierr);
3493 ierr = PCShellSetApply(pc,__mfem_pc_shell_apply); CHKERRQ(ierr);
3494 ierr = PCShellSetApplyTranspose(pc,__mfem_pc_shell_apply_transpose);
3496 ierr = PCShellSetSetUp(pc,__mfem_pc_shell_setup); CHKERRQ(ierr);
3497 ierr = PCShellSetView(pc,__mfem_pc_shell_view); CHKERRQ(ierr);
3498 ierr = PCShellSetDestroy(pc,__mfem_pc_shell_destroy); CHKERRQ(ierr);
3499 PetscFunctionReturn(0);
3504 static PetscErrorCode Convert_Array_IS(MPI_Comm comm,
bool islist,
3506 PetscInt st, IS* is)
3508 PetscInt n = list->
Size(),*idxs;
3509 const int *data = list->
GetData();
3510 PetscErrorCode ierr;
3512 PetscFunctionBeginUser;
3513 ierr = PetscMalloc1(n,&idxs); CHKERRQ(ierr);
3516 for (PetscInt i=0; i<n; i++) { idxs[i] = data[i] + st; }
3521 for (PetscInt i=0; i<n; i++)
3523 if (data[i]) { idxs[cum++] = i+st; }
3527 ierr = ISCreateGeneral(comm,n,idxs,PETSC_OWN_POINTER,is);
3529 PetscFunctionReturn(0);
3535 static PetscErrorCode Convert_Vmarks_IS(MPI_Comm comm,
3538 PetscInt st, IS* is)
3543 PetscErrorCode ierr;
3545 PetscFunctionBeginUser;
3546 for (
int i = 0; i < pl2l.
Size(); i++)
3548 PetscInt m,n,*ii,*jj;
3550 ierr = MatGetRowIJ(pl2l[i],0,PETSC_FALSE,PETSC_FALSE,&m,(
const PetscInt**)&ii,
3551 (
const PetscInt**)&jj,&done); CHKERRQ(ierr);
3552 MFEM_VERIFY(done,
"Unable to perform MatGetRowIJ on " << i <<
" l2l matrix");
3553 ierr = MatGetSize(pl2l[i],NULL,&n); CHKERRQ(ierr);
3557 for (
int i = 0; i < l2l.
Size(); i++) { nl += l2l[i]->Width(); }
3559 const int* vdata = mark->
GetData();
3560 int* sdata = sub_dof_marker.
GetData();
3561 int cumh = 0, cumw = 0;
3562 for (
int i = 0; i < l2l.
Size(); i++)
3567 l2l[i]->BooleanMultTranspose(vf_marker,sf_marker);
3568 cumh += l2l[i]->Height();
3569 cumw += l2l[i]->Width();
3571 ierr = Convert_Array_IS(comm,
false,&sub_dof_marker,st,is); CCHKERRQ(comm,ierr);
3572 for (
int i = 0; i < pl2l.
Size(); i++)
3574 PetscInt m = l2l[i]->Height();
3575 PetscInt *ii = l2l[i]->GetI(),*jj = l2l[i]->GetJ();
3577 ierr = MatRestoreRowIJ(pl2l[i],0,PETSC_FALSE,PETSC_FALSE,&m,
3578 (
const PetscInt**)&ii,
3579 (
const PetscInt**)&jj,&done); CHKERRQ(ierr);
3580 MFEM_VERIFY(done,
"Unable to perform MatRestoreRowIJ on "
3581 << i <<
" l2l matrix");
3584 PetscFunctionReturn(0);
3587 #if !defined(PETSC_HAVE_HYPRE)
3589 #include "_hypre_parcsr_mv.h"
3590 static PetscErrorCode MatConvert_hypreParCSR_AIJ(hypre_ParCSRMatrix* hA,Mat* pA)
3592 MPI_Comm comm = hypre_ParCSRMatrixComm(hA);
3593 hypre_CSRMatrix *hdiag,*hoffd;
3594 PetscScalar *da,*oa,*aptr;
3595 PetscInt *dii,*djj,*oii,*ojj,*iptr;
3596 PetscInt i,dnnz,onnz,m,n;
3598 PetscErrorCode ierr;
3600 PetscFunctionBeginUser;
3601 hdiag = hypre_ParCSRMatrixDiag(hA);
3602 hoffd = hypre_ParCSRMatrixOffd(hA);
3603 m = hypre_CSRMatrixNumRows(hdiag);
3604 n = hypre_CSRMatrixNumCols(hdiag);
3605 dnnz = hypre_CSRMatrixNumNonzeros(hdiag);
3606 onnz = hypre_CSRMatrixNumNonzeros(hoffd);
3607 ierr = PetscMalloc1(m+1,&dii); CHKERRQ(ierr);
3608 ierr = PetscMalloc1(dnnz,&djj); CHKERRQ(ierr);
3609 ierr = PetscMalloc1(dnnz,&da); CHKERRQ(ierr);
3610 ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*
sizeof(PetscInt));
3612 ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*
sizeof(PetscInt));
3614 ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*
sizeof(PetscScalar));
3620 PetscInt nc = dii[i+1]-dii[i];
3621 ierr = PetscSortIntWithScalarArray(nc,iptr,aptr); CHKERRQ(ierr);
3625 ierr = MPI_Comm_size(comm,&size); CHKERRQ(ierr);
3628 PetscInt *offdj,*coffd;
3630 ierr = PetscMalloc1(m+1,&oii); CHKERRQ(ierr);
3631 ierr = PetscMalloc1(onnz,&ojj); CHKERRQ(ierr);
3632 ierr = PetscMalloc1(onnz,&oa); CHKERRQ(ierr);
3633 ierr = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*
sizeof(PetscInt));
3635 offdj = hypre_CSRMatrixJ(hoffd);
3636 coffd = hypre_ParCSRMatrixColMapOffd(hA);
3637 for (i=0; i<onnz; i++) { ojj[i] = coffd[offdj[i]]; }
3638 ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*
sizeof(PetscScalar));
3644 PetscInt nc = oii[i+1]-oii[i];
3645 ierr = PetscSortIntWithScalarArray(nc,iptr,aptr); CHKERRQ(ierr);
3649 ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,
3650 djj,da,oii,ojj,oa,pA); CHKERRQ(ierr);
3656 ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,pA); CHKERRQ(ierr);
3662 void *ptrs[6] = {dii,djj,da,oii,ojj,oa};
3663 const char *names[6] = {
"_mfem_csr_dii",
3674 ierr = PetscContainerCreate(comm,&c); CHKERRQ(ierr);
3675 ierr = PetscContainerSetPointer(c,ptrs[i]); CHKERRQ(ierr);
3676 ierr = PetscContainerSetUserDestroy(c,__mfem_array_container_destroy);
3678 ierr = PetscObjectCompose((PetscObject)(*pA),names[i],(PetscObject)c);
3680 ierr = PetscContainerDestroy(&c); CHKERRQ(ierr);
3682 PetscFunctionReturn(0);
3685 static PetscErrorCode MatConvert_hypreParCSR_IS(hypre_ParCSRMatrix* hA,Mat* pA)
3688 ISLocalToGlobalMapping rl2g,cl2g;
3690 hypre_CSRMatrix *hdiag,*hoffd;
3691 MPI_Comm comm = hypre_ParCSRMatrixComm(hA);
3693 const char *names[2] = {
"_mfem_csr_aux",
3696 PetscScalar *hdd,*hod,*aa,*data;
3697 PetscInt *col_map_offd,*hdi,*hdj,*hoi,*hoj;
3698 PetscInt *aux,*ii,*jj;
3699 PetscInt cum,dr,dc,oc,str,stc,nnz,i,jd,jo;
3700 PetscErrorCode ierr;
3702 PetscFunctionBeginUser;
3704 str = hypre_ParCSRMatrixFirstRowIndex(hA);
3705 stc = hypre_ParCSRMatrixFirstColDiag(hA);
3706 hdiag = hypre_ParCSRMatrixDiag(hA);
3707 hoffd = hypre_ParCSRMatrixOffd(hA);
3708 dr = hypre_CSRMatrixNumRows(hdiag);
3709 dc = hypre_CSRMatrixNumCols(hdiag);
3710 nnz = hypre_CSRMatrixNumNonzeros(hdiag);
3711 hdi = hypre_CSRMatrixI(hdiag);
3712 hdj = hypre_CSRMatrixJ(hdiag);
3713 hdd = hypre_CSRMatrixData(hdiag);
3714 oc = hypre_CSRMatrixNumCols(hoffd);
3715 nnz += hypre_CSRMatrixNumNonzeros(hoffd);
3716 hoi = hypre_CSRMatrixI(hoffd);
3717 hoj = hypre_CSRMatrixJ(hoffd);
3718 hod = hypre_CSRMatrixData(hoffd);
3721 ierr = ISCreateStride(comm,dr,str,1,&is); CHKERRQ(ierr);
3722 ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g); CHKERRQ(ierr);
3723 ierr = ISDestroy(&is); CHKERRQ(ierr);
3724 col_map_offd = hypre_ParCSRMatrixColMapOffd(hA);
3725 ierr = PetscMalloc1(dc+oc,&aux); CHKERRQ(ierr);
3726 for (i=0; i<dc; i++) { aux[i] = i+stc; }
3727 for (i=0; i<oc; i++) { aux[i+dc] = col_map_offd[i]; }
3728 ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is); CHKERRQ(ierr);
3729 ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g); CHKERRQ(ierr);
3730 ierr = ISDestroy(&is); CHKERRQ(ierr);
3733 ierr = MatCreate(comm,pA); CHKERRQ(ierr);
3734 ierr = MatSetSizes(*pA,dr,dc,PETSC_DECIDE,PETSC_DECIDE); CHKERRQ(ierr);
3735 ierr = MatSetType(*pA,MATIS); CHKERRQ(ierr);
3736 ierr = MatSetLocalToGlobalMapping(*pA,rl2g,cl2g); CHKERRQ(ierr);
3737 ierr = ISLocalToGlobalMappingDestroy(&rl2g); CHKERRQ(ierr);
3738 ierr = ISLocalToGlobalMappingDestroy(&cl2g); CHKERRQ(ierr);
3741 ierr = PetscMalloc1(nnz+dr+1,&aux); CHKERRQ(ierr);
3742 ierr = PetscMalloc1(nnz,&data); CHKERRQ(ierr);
3746 *ii = *(hdi++) + *(hoi++);
3747 for (jd=0,jo=0,cum=0; *ii<nnz; cum++)
3749 PetscScalar *aold = aa;
3750 PetscInt *jold = jj,nc = jd+jo;
3751 for (; jd<*hdi; jd++) { *jj++ = *hdj++; *aa++ = *hdd++; }
3752 for (; jo<*hoi; jo++) { *jj++ = *hoj++ + dc; *aa++ = *hod++; }
3753 *(++ii) = *(hdi++) + *(hoi++);
3754 ierr = PetscSortIntWithScalarArray(jd+jo-nc,jold,aold); CHKERRQ(ierr);
3756 for (; cum<dr; cum++) { *(++ii) = nnz; }
3760 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,ii,jj,aa,&lA);
3768 ierr = PetscContainerCreate(PETSC_COMM_SELF,&c); CHKERRQ(ierr);
3769 ierr = PetscContainerSetPointer(c,ptrs[i]); CHKERRQ(ierr);
3770 ierr = PetscContainerSetUserDestroy(c,__mfem_array_container_destroy);
3772 ierr = PetscObjectCompose((PetscObject)lA,names[i],(PetscObject)c);
3774 ierr = PetscContainerDestroy(&c); CHKERRQ(ierr);
3776 ierr = MatISSetLocalMat(*pA,lA); CHKERRQ(ierr);
3777 ierr = MatDestroy(&lA); CHKERRQ(ierr);
3778 ierr = MatAssemblyBegin(*pA,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
3779 ierr = MatAssemblyEnd(*pA,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
3780 PetscFunctionReturn(0);
3784 #endif // MFEM_USE_PETSC
3785 #endif // MFEM_USE_MPI
int Size() const
Logical size of the array.
void SetPreconditioner(Solver &precond)
Sets the solver to perform preconditioning.
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.
Abstract class for PETSc's preconditioners.
PetscParMatrix & operator+=(const PetscParMatrix &B)
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
bool clcustom
Boolean to handle SetFromOptions calls.
void ConvertOperator(MPI_Comm comm, const Operator &op, Mat *B, Operator::Type tid)
Wrapper for PETSc's matrix class.
const Array< int > * nat_dof
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.
virtual void Run(Vector &x, double &t, double &dt, double t_final)
Perform time integration from time t [in] to time tf [in].
ParFiniteElementSpace * fespace
Base abstract class for 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().
Pointer to an Operator of a specified type.
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())
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.
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)
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.
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.
PetscPCGSolver(MPI_Comm comm, const std::string &prefix=std::string())
void Destroy()
Delete all owned data. Does not perform re-initialization with defaults.
virtual Operator & GetExplicitGradient(const Vector &x) const
Return an Operator representing dG/dx at the given point x and the currently set time.
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.
void SortColumnIndices()
Sort the column indices corresponding to each row.
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 Mult(const Vector &b, Vector &x) const
Application of the solver.
virtual ~PetscLinearSolver()
PetscInt GetNumRows() const
Returns the local number of rows.
void Randomize(PetscInt seed)
Set random values.
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()...
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...
virtual void Init(TimeDependentOperator &f_, enum PetscODESolver::Type type)
Initialize the ODE solver.
PetscParMatrix * Transpose(bool action=false)
Returns the transpose of the PetscParMatrix.
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)
void EliminateRowsCols(const Array< int > &rows_cols, const PetscParVector &X, PetscParVector &B)
Eliminate rows and columns from the matrix, and rows from the vector B. Modify B with the BC values i...
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.
virtual Solver * NewPreconditioner(const OperatorHandle &oh)=0
void Assign(const T *)
Copy data from a pointer. Size() elements are copied.
PetscLinearSolver(MPI_Comm comm, const std::string &prefix=std::string(), bool wrap=false)
ID for class PetscParMatrix, MATNEST format.
bool isExplicit() const
True if type is EXPLICIT.
void SetTime(double t)
Sets the current time.
ID for class PetscParMatrix, MATSHELL format.
double * GetData() const
Return element data, i.e. array A.
int * GetI() const
Return the array I.
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.
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.
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...
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.
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...
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.
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.
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.
PetscInt GlobalSize() const
Returns the global number of rows.
const Array< int > * ess_dof
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.
const FiniteElementCollection * FEColl() const
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())
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.
int * GetJ() const
Return the array J.
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.
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.
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...