1674ae819SStefano Zampini #include "bddc.h" 2674ae819SStefano Zampini #include "bddcprivate.h" 3674ae819SStefano Zampini 4674ae819SStefano Zampini #undef __FUNCT__ 5674ae819SStefano Zampini #define __FUNCT__ "PCBDDCNullSpaceAssembleCoarse" 6fdc635d7SStefano Zampini PetscErrorCode PCBDDCNullSpaceAssembleCoarse(PC pc, Mat coarse_mat, MatNullSpace* CoarseNullSpace) 7674ae819SStefano Zampini { 8674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 9674ae819SStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 10*98a51de6SStefano Zampini MatNullSpace tempCoarseNullSpace=NULL; 11674ae819SStefano Zampini const Vec *nsp_vecs; 12*98a51de6SStefano Zampini Vec *coarse_nsp_vecs,local_vec,local_primal_vec,wcoarse_vec,wcoarse_rhs; 13674ae819SStefano Zampini PetscInt nsp_size,coarse_nsp_size,i; 14674ae819SStefano Zampini PetscBool nsp_has_cnst; 15674ae819SStefano Zampini PetscReal test_null; 16674ae819SStefano Zampini PetscErrorCode ierr; 17674ae819SStefano Zampini 18674ae819SStefano Zampini PetscFunctionBegin; 19674ae819SStefano Zampini tempCoarseNullSpace = 0; 20674ae819SStefano Zampini coarse_nsp_size = 0; 21674ae819SStefano Zampini coarse_nsp_vecs = 0; 22674ae819SStefano Zampini ierr = MatNullSpaceGetVecs(pcbddc->NullSpace,&nsp_has_cnst,&nsp_size,&nsp_vecs);CHKERRQ(ierr); 23fdc635d7SStefano Zampini if (coarse_mat) { 24674ae819SStefano Zampini ierr = PetscMalloc((nsp_size+1)*sizeof(Vec),&coarse_nsp_vecs);CHKERRQ(ierr); 25674ae819SStefano Zampini for (i=0;i<nsp_size+1;i++) { 26*98a51de6SStefano Zampini ierr = MatGetVecs(coarse_mat,&coarse_nsp_vecs[i],NULL);CHKERRQ(ierr); 27*98a51de6SStefano Zampini } 28*98a51de6SStefano Zampini if (pcbddc->dbg_flag) { 29*98a51de6SStefano Zampini ierr = MatGetVecs(coarse_mat,&wcoarse_vec,&wcoarse_rhs);CHKERRQ(ierr); 30674ae819SStefano Zampini } 31674ae819SStefano Zampini } 32674ae819SStefano Zampini ierr = MatGetVecs(pcbddc->ConstraintMatrix,&local_vec,&local_primal_vec);CHKERRQ(ierr); 33674ae819SStefano Zampini if (nsp_has_cnst) { 34674ae819SStefano Zampini ierr = VecSet(local_vec,1.0);CHKERRQ(ierr); 35674ae819SStefano Zampini ierr = MatMult(pcbddc->ConstraintMatrix,local_vec,local_primal_vec);CHKERRQ(ierr); 3612edc857SStefano Zampini ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3712edc857SStefano Zampini ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 38fdc635d7SStefano Zampini if (coarse_mat) { 39*98a51de6SStefano Zampini PetscScalar *array_out; 40*98a51de6SStefano Zampini const PetscScalar *array_in; 41*98a51de6SStefano Zampini PetscInt lsize; 42674ae819SStefano Zampini if (pcbddc->dbg_flag) { 43*98a51de6SStefano Zampini PetscViewer dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)coarse_mat)); 44*98a51de6SStefano Zampini ierr = MatMult(coarse_mat,wcoarse_vec,wcoarse_rhs);CHKERRQ(ierr); 45*98a51de6SStefano Zampini ierr = VecNorm(wcoarse_rhs,NORM_INFINITY,&test_null);CHKERRQ(ierr); 46*98a51de6SStefano Zampini ierr = PetscViewerASCIIPrintf(dbg_viewer,"Constant coarse null space error % 1.14e\n",test_null);CHKERRQ(ierr); 47*98a51de6SStefano Zampini ierr = PetscViewerFlush(dbg_viewer);CHKERRQ(ierr); 48674ae819SStefano Zampini } 49*98a51de6SStefano Zampini ierr = VecGetLocalSize(pcbddc->coarse_vec,&lsize);CHKERRQ(ierr); 50*98a51de6SStefano Zampini ierr = VecGetArrayRead(pcbddc->coarse_vec,&array_in);CHKERRQ(ierr); 51*98a51de6SStefano Zampini ierr = VecGetArray(coarse_nsp_vecs[coarse_nsp_size],&array_out);CHKERRQ(ierr); 52*98a51de6SStefano Zampini ierr = PetscMemcpy(array_out,array_in,lsize*sizeof(PetscScalar));CHKERRQ(ierr); 53*98a51de6SStefano Zampini ierr = VecRestoreArray(coarse_nsp_vecs[coarse_nsp_size],&array_out);CHKERRQ(ierr); 54*98a51de6SStefano Zampini ierr = VecRestoreArrayRead(pcbddc->coarse_vec,&array_in);CHKERRQ(ierr); 55674ae819SStefano Zampini coarse_nsp_size++; 56674ae819SStefano Zampini } 57674ae819SStefano Zampini } 58674ae819SStefano Zampini for (i=0;i<nsp_size;i++) { 59674ae819SStefano Zampini ierr = VecScatterBegin(matis->ctx,nsp_vecs[i],local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 60674ae819SStefano Zampini ierr = VecScatterEnd(matis->ctx,nsp_vecs[i],local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 61674ae819SStefano Zampini ierr = MatMult(pcbddc->ConstraintMatrix,local_vec,local_primal_vec);CHKERRQ(ierr); 6212edc857SStefano Zampini ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6312edc857SStefano Zampini ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 64fdc635d7SStefano Zampini if (coarse_mat) { 65*98a51de6SStefano Zampini PetscScalar *array_out; 66*98a51de6SStefano Zampini const PetscScalar *array_in; 67*98a51de6SStefano Zampini PetscInt lsize; 68674ae819SStefano Zampini if (pcbddc->dbg_flag) { 69*98a51de6SStefano Zampini PetscViewer dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)coarse_mat)); 70*98a51de6SStefano Zampini ierr = MatMult(coarse_mat,wcoarse_vec,wcoarse_rhs);CHKERRQ(ierr); 71*98a51de6SStefano Zampini ierr = VecNorm(wcoarse_rhs,NORM_2,&test_null);CHKERRQ(ierr); 72*98a51de6SStefano Zampini ierr = PetscViewerASCIIPrintf(dbg_viewer,"Vec %d coarse null space error % 1.14e\n",i,test_null);CHKERRQ(ierr); 73*98a51de6SStefano Zampini ierr = PetscViewerFlush(dbg_viewer);CHKERRQ(ierr); 74674ae819SStefano Zampini } 75*98a51de6SStefano Zampini ierr = VecGetLocalSize(pcbddc->coarse_vec,&lsize);CHKERRQ(ierr); 76*98a51de6SStefano Zampini ierr = VecGetArrayRead(pcbddc->coarse_vec,&array_in);CHKERRQ(ierr); 77*98a51de6SStefano Zampini ierr = VecGetArray(coarse_nsp_vecs[coarse_nsp_size],&array_out);CHKERRQ(ierr); 78*98a51de6SStefano Zampini ierr = PetscMemcpy(array_out,array_in,lsize*sizeof(PetscScalar));CHKERRQ(ierr); 79*98a51de6SStefano Zampini ierr = VecRestoreArray(coarse_nsp_vecs[coarse_nsp_size],&array_out);CHKERRQ(ierr); 80*98a51de6SStefano Zampini ierr = VecRestoreArrayRead(pcbddc->coarse_vec,&array_in);CHKERRQ(ierr); 81674ae819SStefano Zampini coarse_nsp_size++; 82674ae819SStefano Zampini } 83674ae819SStefano Zampini } 84674ae819SStefano Zampini if (coarse_nsp_size > 0) { 859a7d3425SStefano Zampini ierr = PCBDDCOrthonormalizeVecs(coarse_nsp_size,coarse_nsp_vecs);CHKERRQ(ierr); 86fdc635d7SStefano Zampini ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject)coarse_mat),PETSC_FALSE,coarse_nsp_size,coarse_nsp_vecs,&tempCoarseNullSpace);CHKERRQ(ierr); 87674ae819SStefano Zampini for (i=0;i<nsp_size+1;i++) { 88674ae819SStefano Zampini ierr = VecDestroy(&coarse_nsp_vecs[i]);CHKERRQ(ierr); 89674ae819SStefano Zampini } 90674ae819SStefano Zampini } 91*98a51de6SStefano Zampini if (coarse_mat) { 92674ae819SStefano Zampini ierr = PetscFree(coarse_nsp_vecs);CHKERRQ(ierr); 93*98a51de6SStefano Zampini if (pcbddc->dbg_flag) { 94*98a51de6SStefano Zampini ierr = VecDestroy(&wcoarse_vec);CHKERRQ(ierr); 95*98a51de6SStefano Zampini ierr = VecDestroy(&wcoarse_rhs);CHKERRQ(ierr); 96*98a51de6SStefano Zampini } 97*98a51de6SStefano Zampini } 98674ae819SStefano Zampini ierr = VecDestroy(&local_vec);CHKERRQ(ierr); 99674ae819SStefano Zampini ierr = VecDestroy(&local_primal_vec);CHKERRQ(ierr); 100674ae819SStefano Zampini *CoarseNullSpace = tempCoarseNullSpace; 101674ae819SStefano Zampini PetscFunctionReturn(0); 102674ae819SStefano Zampini } 103674ae819SStefano Zampini 104674ae819SStefano Zampini 105674ae819SStefano Zampini #undef __FUNCT__ 106674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplyNullSpaceCorrectionPC" 107674ae819SStefano Zampini static PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC pc,Vec x,Vec y) 108674ae819SStefano Zampini { 109674ae819SStefano Zampini NullSpaceCorrection_ctx pc_ctx; 110674ae819SStefano Zampini PetscErrorCode ierr; 111674ae819SStefano Zampini 112674ae819SStefano Zampini PetscFunctionBegin; 113674ae819SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 114674ae819SStefano Zampini /* E */ 115674ae819SStefano Zampini ierr = MatMultTranspose(pc_ctx->Lbasis_mat,x,pc_ctx->work_small_2);CHKERRQ(ierr); 116674ae819SStefano Zampini ierr = MatMultAdd(pc_ctx->Kbasis_mat,pc_ctx->work_small_2,x,pc_ctx->work_full_1);CHKERRQ(ierr); 117674ae819SStefano Zampini /* P^-1 */ 118674ae819SStefano Zampini ierr = PCApply(pc_ctx->local_pc,pc_ctx->work_full_1,pc_ctx->work_full_2);CHKERRQ(ierr); 119674ae819SStefano Zampini /* E^T */ 120674ae819SStefano Zampini ierr = MatMultTranspose(pc_ctx->Kbasis_mat,pc_ctx->work_full_2,pc_ctx->work_small_1);CHKERRQ(ierr); 121674ae819SStefano Zampini ierr = VecScale(pc_ctx->work_small_1,-1.0);CHKERRQ(ierr); 122674ae819SStefano Zampini ierr = MatMultAdd(pc_ctx->Lbasis_mat,pc_ctx->work_small_1,pc_ctx->work_full_2,pc_ctx->work_full_1);CHKERRQ(ierr); 123674ae819SStefano Zampini /* Sum contributions */ 124674ae819SStefano Zampini ierr = MatMultAdd(pc_ctx->basis_mat,pc_ctx->work_small_2,pc_ctx->work_full_1,y);CHKERRQ(ierr); 125674ae819SStefano Zampini PetscFunctionReturn(0); 126674ae819SStefano Zampini } 127674ae819SStefano Zampini 128674ae819SStefano Zampini #undef __FUNCT__ 129674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyNullSpaceCorrectionPC" 130674ae819SStefano Zampini static PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC pc) 131674ae819SStefano Zampini { 132674ae819SStefano Zampini NullSpaceCorrection_ctx pc_ctx; 133674ae819SStefano Zampini PetscErrorCode ierr; 134674ae819SStefano Zampini 135674ae819SStefano Zampini PetscFunctionBegin; 136674ae819SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 137674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_small_1);CHKERRQ(ierr); 138674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_small_2);CHKERRQ(ierr); 139674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_full_1);CHKERRQ(ierr); 140674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_full_2);CHKERRQ(ierr); 141674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->basis_mat);CHKERRQ(ierr); 142674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->Lbasis_mat);CHKERRQ(ierr); 143674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->Kbasis_mat);CHKERRQ(ierr); 144674ae819SStefano Zampini ierr = PCDestroy(&pc_ctx->local_pc);CHKERRQ(ierr); 145674ae819SStefano Zampini ierr = PetscFree(pc_ctx);CHKERRQ(ierr); 146674ae819SStefano Zampini PetscFunctionReturn(0); 147674ae819SStefano Zampini } 148674ae819SStefano Zampini 14928509bceSStefano Zampini /* 15028509bceSStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC,Vec,Vec); 15128509bceSStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC); 15228509bceSStefano Zampini */ 153674ae819SStefano Zampini 154674ae819SStefano Zampini #undef __FUNCT__ 155674ae819SStefano Zampini #define __FUNCT__ "PCBDDCNullSpaceAssembleCorrection" 156674ae819SStefano Zampini PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc,IS local_dofs) 157674ae819SStefano Zampini { 158674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 159674ae819SStefano Zampini PC_IS *pcis = (PC_IS*)pc->data; 160674ae819SStefano Zampini Mat_IS* matis = (Mat_IS*)pc->pmat->data; 161674ae819SStefano Zampini KSP *local_ksp; 162674ae819SStefano Zampini PC newpc; 163674ae819SStefano Zampini NullSpaceCorrection_ctx shell_ctx; 164674ae819SStefano Zampini Mat local_mat,local_pmat,small_mat,inv_small_mat; 165674ae819SStefano Zampini MatStructure local_mat_struct; 166674ae819SStefano Zampini Vec work1,work2; 167674ae819SStefano Zampini const Vec *nullvecs; 168674ae819SStefano Zampini VecScatter scatter_ctx; 169674ae819SStefano Zampini IS is_aux; 170674ae819SStefano Zampini MatFactorInfo matinfo; 171674ae819SStefano Zampini PetscScalar *basis_mat,*Kbasis_mat,*array,*array_mat; 172674ae819SStefano Zampini PetscScalar one = 1.0,zero = 0.0, m_one = -1.0; 173674ae819SStefano Zampini PetscInt basis_dofs,basis_size,nnsp_size,i,k,n_I,n_R; 174674ae819SStefano Zampini PetscBool nnsp_has_cnst; 175674ae819SStefano Zampini PetscErrorCode ierr; 176674ae819SStefano Zampini 177674ae819SStefano Zampini PetscFunctionBegin; 178674ae819SStefano Zampini /* Infer the local solver */ 179674ae819SStefano Zampini ierr = ISGetSize(local_dofs,&basis_dofs);CHKERRQ(ierr); 180674ae819SStefano Zampini ierr = VecGetSize(pcis->vec1_D,&n_I);CHKERRQ(ierr); 181674ae819SStefano Zampini ierr = VecGetSize(pcbddc->vec1_R,&n_R);CHKERRQ(ierr); 182674ae819SStefano Zampini if (basis_dofs == n_I) { 183674ae819SStefano Zampini /* Dirichlet solver */ 184674ae819SStefano Zampini local_ksp = &pcbddc->ksp_D; 185674ae819SStefano Zampini } else if (basis_dofs == n_R) { 186674ae819SStefano Zampini /* Neumann solver */ 187674ae819SStefano Zampini local_ksp = &pcbddc->ksp_R; 188674ae819SStefano Zampini } else { 189674ae819SStefano Zampini SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in %s: unknown local IS size %d. n_I=%d, n_R=%d)\n",__FUNCT__,basis_dofs,n_I,n_R); 190674ae819SStefano Zampini } 191674ae819SStefano Zampini ierr = KSPGetOperators(*local_ksp,&local_mat,&local_pmat,&local_mat_struct);CHKERRQ(ierr); 192674ae819SStefano Zampini 193674ae819SStefano Zampini /* Get null space vecs */ 194674ae819SStefano Zampini ierr = MatNullSpaceGetVecs(pcbddc->NullSpace,&nnsp_has_cnst,&nnsp_size,&nullvecs);CHKERRQ(ierr); 195674ae819SStefano Zampini basis_size = nnsp_size; 196674ae819SStefano Zampini if (nnsp_has_cnst) { 197674ae819SStefano Zampini basis_size++; 198674ae819SStefano Zampini } 199674ae819SStefano Zampini 200b8ffe317SStefano Zampini if (basis_dofs) { 201674ae819SStefano Zampini /* Create shell ctx */ 202674ae819SStefano Zampini ierr = PetscMalloc(sizeof(*shell_ctx),&shell_ctx);CHKERRQ(ierr); 203674ae819SStefano Zampini 204674ae819SStefano Zampini /* Create work vectors in shell context */ 205674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_small_1);CHKERRQ(ierr); 206674ae819SStefano Zampini ierr = VecSetSizes(shell_ctx->work_small_1,basis_size,basis_size);CHKERRQ(ierr); 207674ae819SStefano Zampini ierr = VecSetType(shell_ctx->work_small_1,VECSEQ);CHKERRQ(ierr); 208674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_small_1,&shell_ctx->work_small_2);CHKERRQ(ierr); 209674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_full_1);CHKERRQ(ierr); 210674ae819SStefano Zampini ierr = VecSetSizes(shell_ctx->work_full_1,basis_dofs,basis_dofs);CHKERRQ(ierr); 211674ae819SStefano Zampini ierr = VecSetType(shell_ctx->work_full_1,VECSEQ);CHKERRQ(ierr); 212674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&shell_ctx->work_full_2);CHKERRQ(ierr); 213674ae819SStefano Zampini 214674ae819SStefano Zampini /* Allocate workspace */ 215674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->basis_mat );CHKERRQ(ierr); 216674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->Kbasis_mat);CHKERRQ(ierr); 217674ae819SStefano Zampini ierr = MatDenseGetArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr); 218674ae819SStefano Zampini ierr = MatDenseGetArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr); 219674ae819SStefano Zampini 220674ae819SStefano Zampini /* Restrict local null space on selected dofs (Dirichlet or Neumann) 221674ae819SStefano Zampini and compute matrices N and K*N */ 222674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work1);CHKERRQ(ierr); 223674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work2);CHKERRQ(ierr); 224674ae819SStefano Zampini ierr = VecScatterCreate(pcis->vec1_N,local_dofs,work1,(IS)0,&scatter_ctx);CHKERRQ(ierr); 225b8ffe317SStefano Zampini } 226b8ffe317SStefano Zampini 227674ae819SStefano Zampini for (k=0;k<nnsp_size;k++) { 228674ae819SStefano Zampini ierr = VecScatterBegin(matis->ctx,nullvecs[k],pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 229674ae819SStefano Zampini ierr = VecScatterEnd(matis->ctx,nullvecs[k],pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 230b8ffe317SStefano Zampini if (basis_dofs) { 231674ae819SStefano Zampini ierr = VecPlaceArray(work1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr); 232674ae819SStefano Zampini ierr = VecScatterBegin(scatter_ctx,pcis->vec1_N,work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 233674ae819SStefano Zampini ierr = VecScatterEnd(scatter_ctx,pcis->vec1_N,work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 234674ae819SStefano Zampini ierr = VecPlaceArray(work2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr); 235674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr); 236674ae819SStefano Zampini ierr = VecResetArray(work1);CHKERRQ(ierr); 237674ae819SStefano Zampini ierr = VecResetArray(work2);CHKERRQ(ierr); 238674ae819SStefano Zampini } 239b8ffe317SStefano Zampini } 240b8ffe317SStefano Zampini 241b8ffe317SStefano Zampini if (basis_dofs) { 242674ae819SStefano Zampini if (nnsp_has_cnst) { 243674ae819SStefano Zampini ierr = VecPlaceArray(work1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr); 244674ae819SStefano Zampini ierr = VecSet(work1,one);CHKERRQ(ierr); 245674ae819SStefano Zampini ierr = VecPlaceArray(work2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr); 246674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr); 247674ae819SStefano Zampini ierr = VecResetArray(work1);CHKERRQ(ierr); 248674ae819SStefano Zampini ierr = VecResetArray(work2);CHKERRQ(ierr); 249674ae819SStefano Zampini } 250674ae819SStefano Zampini ierr = VecDestroy(&work1);CHKERRQ(ierr); 251674ae819SStefano Zampini ierr = VecDestroy(&work2);CHKERRQ(ierr); 252674ae819SStefano Zampini ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr); 253674ae819SStefano Zampini ierr = MatDenseRestoreArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr); 254674ae819SStefano Zampini ierr = MatDenseRestoreArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr); 255674ae819SStefano Zampini 256674ae819SStefano Zampini /* Assemble another Mat object in shell context */ 257674ae819SStefano Zampini ierr = MatTransposeMatMult(shell_ctx->basis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&small_mat);CHKERRQ(ierr); 258674ae819SStefano Zampini ierr = MatFactorInfoInitialize(&matinfo);CHKERRQ(ierr); 259674ae819SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,basis_size,0,1,&is_aux);CHKERRQ(ierr); 260674ae819SStefano Zampini ierr = MatLUFactor(small_mat,is_aux,is_aux,&matinfo);CHKERRQ(ierr); 261674ae819SStefano Zampini ierr = ISDestroy(&is_aux);CHKERRQ(ierr); 262674ae819SStefano Zampini ierr = PetscMalloc(basis_size*basis_size*sizeof(PetscScalar),&array_mat);CHKERRQ(ierr); 263674ae819SStefano Zampini for (k=0;k<basis_size;k++) { 264674ae819SStefano Zampini ierr = VecSet(shell_ctx->work_small_1,zero);CHKERRQ(ierr); 265674ae819SStefano Zampini ierr = VecSetValue(shell_ctx->work_small_1,k,one,INSERT_VALUES);CHKERRQ(ierr); 266674ae819SStefano Zampini ierr = VecAssemblyBegin(shell_ctx->work_small_1);CHKERRQ(ierr); 267674ae819SStefano Zampini ierr = VecAssemblyEnd(shell_ctx->work_small_1);CHKERRQ(ierr); 268674ae819SStefano Zampini ierr = MatSolve(small_mat,shell_ctx->work_small_1,shell_ctx->work_small_2);CHKERRQ(ierr); 269674ae819SStefano Zampini ierr = VecGetArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr); 270674ae819SStefano Zampini for (i=0;i<basis_size;i++) { 271674ae819SStefano Zampini array_mat[i*basis_size+k]=array[i]; 272674ae819SStefano Zampini } 273674ae819SStefano Zampini ierr = VecRestoreArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr); 274674ae819SStefano Zampini } 275674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_size,basis_size,array_mat,&inv_small_mat);CHKERRQ(ierr); 276674ae819SStefano Zampini ierr = MatMatMult(shell_ctx->basis_mat,inv_small_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->Lbasis_mat);CHKERRQ(ierr); 277674ae819SStefano Zampini ierr = PetscFree(array_mat);CHKERRQ(ierr); 278674ae819SStefano Zampini ierr = MatDestroy(&inv_small_mat);CHKERRQ(ierr); 279674ae819SStefano Zampini ierr = MatDestroy(&small_mat);CHKERRQ(ierr); 280674ae819SStefano Zampini ierr = MatScale(shell_ctx->Kbasis_mat,m_one);CHKERRQ(ierr); 281674ae819SStefano Zampini 282674ae819SStefano Zampini /* Rebuild local PC */ 283674ae819SStefano Zampini ierr = KSPGetPC(*local_ksp,&shell_ctx->local_pc);CHKERRQ(ierr); 284674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)shell_ctx->local_pc);CHKERRQ(ierr); 285674ae819SStefano Zampini ierr = PCCreate(PETSC_COMM_SELF,&newpc);CHKERRQ(ierr); 286674ae819SStefano Zampini ierr = PCSetOperators(newpc,local_mat,local_mat,SAME_PRECONDITIONER);CHKERRQ(ierr); 287674ae819SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 288674ae819SStefano Zampini ierr = PCShellSetContext(newpc,shell_ctx);CHKERRQ(ierr); 289674ae819SStefano Zampini ierr = PCShellSetApply(newpc,PCBDDCApplyNullSpaceCorrectionPC);CHKERRQ(ierr); 290674ae819SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyNullSpaceCorrectionPC);CHKERRQ(ierr); 291674ae819SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 292674ae819SStefano Zampini ierr = KSPSetPC(*local_ksp,newpc);CHKERRQ(ierr); 293674ae819SStefano Zampini ierr = PCDestroy(&newpc);CHKERRQ(ierr); 294674ae819SStefano Zampini ierr = KSPSetUp(*local_ksp);CHKERRQ(ierr); 295b8ffe317SStefano Zampini } 296674ae819SStefano Zampini /* test */ 297b8ffe317SStefano Zampini if (pcbddc->dbg_flag && basis_dofs) { 298674ae819SStefano Zampini KSP check_ksp; 299674ae819SStefano Zampini PC check_pc; 300674ae819SStefano Zampini Mat test_mat; 301674ae819SStefano Zampini Vec work3; 302674ae819SStefano Zampini PetscReal test_err,lambda_min,lambda_max; 303674ae819SStefano Zampini PetscBool setsym,issym=PETSC_FALSE; 304b8ffe317SStefano Zampini PetscInt tabs; 305674ae819SStefano Zampini 306b8ffe317SStefano Zampini ierr = PetscViewerASCIIGetTab(pcbddc->dbg_viewer,&tabs);CHKERRQ(ierr); 307674ae819SStefano Zampini ierr = KSPGetPC(*local_ksp,&check_pc);CHKERRQ(ierr); 308674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work1);CHKERRQ(ierr); 309674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work2);CHKERRQ(ierr); 310674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work3);CHKERRQ(ierr); 311674ae819SStefano Zampini ierr = VecSetRandom(shell_ctx->work_small_1,NULL);CHKERRQ(ierr); 312674ae819SStefano Zampini ierr = MatMult(shell_ctx->basis_mat,shell_ctx->work_small_1,work1);CHKERRQ(ierr); 313674ae819SStefano Zampini ierr = VecCopy(work1,work2);CHKERRQ(ierr); 314674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work3);CHKERRQ(ierr); 315674ae819SStefano Zampini ierr = PCApply(check_pc,work3,work1);CHKERRQ(ierr); 316674ae819SStefano Zampini ierr = VecAXPY(work1,m_one,work2);CHKERRQ(ierr); 317674ae819SStefano Zampini ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr); 318b8ffe317SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d error for nullspace correction for ",PetscGlobalRank); 319b8ffe317SStefano Zampini ierr = PetscViewerASCIIUseTabs(pcbddc->dbg_viewer,PETSC_FALSE);CHKERRQ(ierr); 320674ae819SStefano Zampini if (basis_dofs == n_I) { 321b8ffe317SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Dirichlet "); 322674ae819SStefano Zampini } else { 323b8ffe317SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Neumann "); 324674ae819SStefano Zampini } 325b8ffe317SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"solver is :%1.14e\n",test_err); 326b8ffe317SStefano Zampini ierr = PetscViewerASCIISetTab(pcbddc->dbg_viewer,tabs);CHKERRQ(ierr); 327b8ffe317SStefano Zampini ierr = PetscViewerASCIIUseTabs(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 328674ae819SStefano Zampini 329674ae819SStefano Zampini ierr = MatTransposeMatMult(shell_ctx->Lbasis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&test_mat);CHKERRQ(ierr); 330674ae819SStefano Zampini ierr = MatShift(test_mat,one);CHKERRQ(ierr); 331674ae819SStefano Zampini ierr = MatNorm(test_mat,NORM_INFINITY,&test_err);CHKERRQ(ierr); 332674ae819SStefano Zampini ierr = MatDestroy(&test_mat);CHKERRQ(ierr); 333b8ffe317SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d error for nullspace matrices is :%1.14e\n",PetscGlobalRank,test_err); 334674ae819SStefano Zampini 335674ae819SStefano Zampini /* Create ksp object suitable for extreme eigenvalues' estimation */ 336674ae819SStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr); 337674ae819SStefano Zampini ierr = KSPSetOperators(check_ksp,local_mat,local_mat,SAME_PRECONDITIONER);CHKERRQ(ierr); 338674ae819SStefano Zampini ierr = KSPSetTolerances(check_ksp,1.e-8,1.e-8,PETSC_DEFAULT,basis_dofs);CHKERRQ(ierr); 339674ae819SStefano Zampini ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr); 340674ae819SStefano Zampini ierr = MatIsSymmetricKnown(pc->pmat,&setsym,&issym);CHKERRQ(ierr); 341674ae819SStefano Zampini if (issym) { 342674ae819SStefano Zampini ierr = KSPSetType(check_ksp,KSPCG);CHKERRQ(ierr); 343674ae819SStefano Zampini } 344674ae819SStefano Zampini ierr = KSPSetPC(check_ksp,check_pc);CHKERRQ(ierr); 345674ae819SStefano Zampini ierr = KSPSetUp(check_ksp);CHKERRQ(ierr); 346674ae819SStefano Zampini ierr = VecSetRandom(work1,NULL);CHKERRQ(ierr); 347674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr); 348674ae819SStefano Zampini ierr = KSPSolve(check_ksp,work2,work2);CHKERRQ(ierr); 349674ae819SStefano Zampini ierr = VecAXPY(work2,m_one,work1);CHKERRQ(ierr); 350674ae819SStefano Zampini ierr = VecNorm(work2,NORM_INFINITY,&test_err);CHKERRQ(ierr); 351674ae819SStefano Zampini ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr); 352674ae819SStefano Zampini ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr); 353b8ffe317SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d error for adapted KSP %1.14e (it %d, eigs %1.6e %1.6e)\n",PetscGlobalRank,test_err,k,lambda_min,lambda_max); 354674ae819SStefano Zampini ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr); 355674ae819SStefano Zampini ierr = VecDestroy(&work1);CHKERRQ(ierr); 356674ae819SStefano Zampini ierr = VecDestroy(&work2);CHKERRQ(ierr); 357674ae819SStefano Zampini ierr = VecDestroy(&work3);CHKERRQ(ierr); 358674ae819SStefano Zampini } 359b8ffe317SStefano Zampini /* all processes shoud call this, even the void ones */ 360b8ffe317SStefano Zampini if (pcbddc->dbg_flag) { 361b8ffe317SStefano Zampini ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 362b8ffe317SStefano Zampini } 363674ae819SStefano Zampini PetscFunctionReturn(0); 364674ae819SStefano Zampini } 365674ae819SStefano Zampini 366*98a51de6SStefano Zampini /* uncomment to test nullspace adaptation when change of basis has been requested */ 367*98a51de6SStefano Zampini /* #define PCBDDC_TESTNULLSPACE 1 */ 368674ae819SStefano Zampini #undef __FUNCT__ 369674ae819SStefano Zampini #define __FUNCT__ "PCBDDCNullSpaceAdaptGlobal" 370674ae819SStefano Zampini PetscErrorCode PCBDDCNullSpaceAdaptGlobal(PC pc) 371674ae819SStefano Zampini { 372674ae819SStefano Zampini PC_IS* pcis = (PC_IS*)(pc->data); 373674ae819SStefano Zampini PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 374674ae819SStefano Zampini KSP inv_change; 375674ae819SStefano Zampini PC pc_change; 376674ae819SStefano Zampini const Vec *nsp_vecs; 377674ae819SStefano Zampini Vec *new_nsp_vecs; 378674ae819SStefano Zampini PetscInt i,nsp_size,new_nsp_size,start_new; 379674ae819SStefano Zampini PetscBool nsp_has_cnst; 380674ae819SStefano Zampini MatNullSpace new_nsp; 381674ae819SStefano Zampini PetscErrorCode ierr; 382674ae819SStefano Zampini 383674ae819SStefano Zampini PetscFunctionBegin; 384674ae819SStefano Zampini /* create KSP for change of basis */ 385674ae819SStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&inv_change);CHKERRQ(ierr); 386674ae819SStefano Zampini ierr = KSPSetOperators(inv_change,pcbddc->ChangeOfBasisMatrix,pcbddc->ChangeOfBasisMatrix,SAME_PRECONDITIONER);CHKERRQ(ierr); 387674ae819SStefano Zampini ierr = KSPSetType(inv_change,KSPPREONLY);CHKERRQ(ierr); 388674ae819SStefano Zampini ierr = KSPGetPC(inv_change,&pc_change);CHKERRQ(ierr); 389674ae819SStefano Zampini ierr = PCSetType(pc_change,PCLU);CHKERRQ(ierr); 390674ae819SStefano Zampini ierr = KSPSetUp(inv_change);CHKERRQ(ierr); 391674ae819SStefano Zampini /* get nullspace and transform it */ 392674ae819SStefano Zampini ierr = MatNullSpaceGetVecs(pcbddc->NullSpace,&nsp_has_cnst,&nsp_size,&nsp_vecs);CHKERRQ(ierr); 393674ae819SStefano Zampini new_nsp_size = nsp_size; 394674ae819SStefano Zampini if (nsp_has_cnst) { 395674ae819SStefano Zampini new_nsp_size++; 396674ae819SStefano Zampini } 397674ae819SStefano Zampini ierr = PetscMalloc(new_nsp_size*sizeof(Vec),&new_nsp_vecs);CHKERRQ(ierr); 398674ae819SStefano Zampini for (i=0;i<new_nsp_size;i++) { 399674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&new_nsp_vecs[i]);CHKERRQ(ierr); 400674ae819SStefano Zampini } 401674ae819SStefano Zampini start_new = 0; 402674ae819SStefano Zampini if (nsp_has_cnst) { 403674ae819SStefano Zampini start_new = 1; 404674ae819SStefano Zampini ierr = VecSet(new_nsp_vecs[0],1.0);CHKERRQ(ierr); 405674ae819SStefano Zampini ierr = VecSet(pcis->vec1_B,1.0);CHKERRQ(ierr); 406*98a51de6SStefano Zampini ierr = KSPSolve(inv_change,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr); 407674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[0],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 408674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[0],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 409674ae819SStefano Zampini } 410674ae819SStefano Zampini for (i=0;i<nsp_size;i++) { 411674ae819SStefano Zampini ierr = VecCopy(nsp_vecs[i],new_nsp_vecs[i+start_new]);CHKERRQ(ierr); 412674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,nsp_vecs[i],pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 413674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,nsp_vecs[i],pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 414674ae819SStefano Zampini ierr = KSPSolve(inv_change,pcis->vec1_B,pcis->vec1_B); 415674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[i+start_new],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 416674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[i+start_new],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 417674ae819SStefano Zampini } 4189a7d3425SStefano Zampini ierr = PCBDDCOrthonormalizeVecs(new_nsp_size,new_nsp_vecs);CHKERRQ(ierr); 419*98a51de6SStefano Zampini #ifdef PCBDDC_TESTNULLSPACE 420*98a51de6SStefano Zampini { 421674ae819SStefano Zampini PetscBool nsp_t=PETSC_FALSE; 422*98a51de6SStefano Zampini PetscReal test_norm; 423674ae819SStefano Zampini Mat temp_mat; 424674ae819SStefano Zampini Mat_IS* matis = (Mat_IS*)pc->pmat->data; 425*98a51de6SStefano Zampini ierr = PetscPrintf(PetscObjectComm((PetscObject)(pc->pmat)),"Testing BDDC nullspace (mat nullspace)\n",nsp_t); 426*98a51de6SStefano Zampini ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr); 427*98a51de6SStefano Zampini ierr = PetscPrintf(PetscObjectComm((PetscObject)(pc->pmat)),"ORI ORI: %d\n",nsp_t); 428674ae819SStefano Zampini temp_mat = matis->A; 429674ae819SStefano Zampini matis->A = pcbddc->local_mat; 430674ae819SStefano Zampini pcbddc->local_mat = temp_mat; 431674ae819SStefano Zampini ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr); 432*98a51de6SStefano Zampini ierr = PetscPrintf(PetscObjectComm((PetscObject)(pc->pmat)),"NEW ORI: %d\n",nsp_t); 433674ae819SStefano Zampini for (i=0;i<new_nsp_size;i++) { 434674ae819SStefano Zampini ierr = MatMult(pc->pmat,new_nsp_vecs[i],pcis->vec1_global);CHKERRQ(ierr); 435674ae819SStefano Zampini ierr = VecNorm(pcis->vec1_global,NORM_2,&test_norm);CHKERRQ(ierr); 436674ae819SStefano Zampini if (test_norm > 1.e-12) { 437*98a51de6SStefano Zampini ierr = PetscPrintf(PetscObjectComm((PetscObject)(pc->pmat)),"------------ERROR VEC %d------------------\n",i); 438674ae819SStefano Zampini ierr = VecView(pcis->vec1_global,PETSC_VIEWER_STDOUT_WORLD); 439*98a51de6SStefano Zampini ierr = PetscPrintf(PetscObjectComm((PetscObject)(pc->pmat)),"------------------------------------------\n"); 440674ae819SStefano Zampini } 441674ae819SStefano Zampini } 442674ae819SStefano Zampini } 443674ae819SStefano Zampini #endif 444674ae819SStefano Zampini 445674ae819SStefano Zampini ierr = KSPDestroy(&inv_change);CHKERRQ(ierr); 446*98a51de6SStefano Zampini ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject)pc),PETSC_FALSE,new_nsp_size,new_nsp_vecs,&new_nsp);CHKERRQ(ierr); 447674ae819SStefano Zampini ierr = PCBDDCSetNullSpace(pc,new_nsp);CHKERRQ(ierr); 448674ae819SStefano Zampini ierr = MatNullSpaceDestroy(&new_nsp);CHKERRQ(ierr); 449*98a51de6SStefano Zampini #ifdef PCBDDC_TESTNULLSPACE 450*98a51de6SStefano Zampini { 451*98a51de6SStefano Zampini PetscBool nsp_t=PETSC_FALSE; 452*98a51de6SStefano Zampini Mat temp_mat; 453*98a51de6SStefano Zampini Mat_IS* matis = (Mat_IS*)pc->pmat->data; 454674ae819SStefano Zampini ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr); 455*98a51de6SStefano Zampini ierr = PetscPrintf(PetscObjectComm((PetscObject)(pc->pmat)),"NEW NEW: %d\n",nsp_t); 456674ae819SStefano Zampini temp_mat = matis->A; 457674ae819SStefano Zampini matis->A = pcbddc->local_mat; 458674ae819SStefano Zampini pcbddc->local_mat = temp_mat; 459674ae819SStefano Zampini ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr); 460*98a51de6SStefano Zampini ierr = PetscPrintf(PetscObjectComm((PetscObject)(pc->pmat)),"ORI NEW: %d\n",nsp_t); 461*98a51de6SStefano Zampini } 462674ae819SStefano Zampini #endif 463674ae819SStefano Zampini 464674ae819SStefano Zampini for (i=0;i<new_nsp_size;i++) { 465674ae819SStefano Zampini ierr = VecDestroy(&new_nsp_vecs[i]);CHKERRQ(ierr); 466674ae819SStefano Zampini } 467674ae819SStefano Zampini ierr = PetscFree(new_nsp_vecs);CHKERRQ(ierr); 468674ae819SStefano Zampini PetscFunctionReturn(0); 469674ae819SStefano Zampini } 470