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; 10674ae819SStefano Zampini MatNullSpace tempCoarseNullSpace; 11674ae819SStefano Zampini const Vec *nsp_vecs; 12674ae819SStefano Zampini Vec *coarse_nsp_vecs,local_vec,local_primal_vec; 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) { 24*785e854fSJed Brown ierr = PetscMalloc1((nsp_size+1),&coarse_nsp_vecs);CHKERRQ(ierr); 25674ae819SStefano Zampini for (i=0;i<nsp_size+1;i++) { 26674ae819SStefano Zampini ierr = VecDuplicate(pcbddc->coarse_vec,&coarse_nsp_vecs[i]);CHKERRQ(ierr); 27674ae819SStefano Zampini } 28674ae819SStefano Zampini } 29674ae819SStefano Zampini ierr = MatGetVecs(pcbddc->ConstraintMatrix,&local_vec,&local_primal_vec);CHKERRQ(ierr); 30674ae819SStefano Zampini if (nsp_has_cnst) { 31674ae819SStefano Zampini ierr = VecSet(local_vec,1.0);CHKERRQ(ierr); 32674ae819SStefano Zampini ierr = MatMult(pcbddc->ConstraintMatrix,local_vec,local_primal_vec);CHKERRQ(ierr); 33674ae819SStefano Zampini ierr = PCBDDCScatterCoarseDataBegin(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 34674ae819SStefano Zampini ierr = PCBDDCScatterCoarseDataEnd(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 35fdc635d7SStefano Zampini if (coarse_mat) { 36674ae819SStefano Zampini if (pcbddc->dbg_flag) { 37fdc635d7SStefano Zampini ierr = MatMult(coarse_mat,pcbddc->coarse_vec,pcbddc->coarse_rhs);CHKERRQ(ierr); 38674ae819SStefano Zampini ierr = VecNorm(pcbddc->coarse_rhs,NORM_INFINITY,&test_null);CHKERRQ(ierr); 39674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Constant coarse null space error % 1.14e\n",test_null);CHKERRQ(ierr); 40674ae819SStefano Zampini } 41674ae819SStefano Zampini ierr = VecCopy(pcbddc->coarse_vec,coarse_nsp_vecs[coarse_nsp_size]);CHKERRQ(ierr); 42674ae819SStefano Zampini coarse_nsp_size++; 43674ae819SStefano Zampini } 44674ae819SStefano Zampini } 45674ae819SStefano Zampini for (i=0;i<nsp_size;i++) { 46674ae819SStefano Zampini ierr = VecScatterBegin(matis->ctx,nsp_vecs[i],local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 47674ae819SStefano Zampini ierr = VecScatterEnd(matis->ctx,nsp_vecs[i],local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 48674ae819SStefano Zampini ierr = MatMult(pcbddc->ConstraintMatrix,local_vec,local_primal_vec);CHKERRQ(ierr); 49674ae819SStefano Zampini ierr = PCBDDCScatterCoarseDataBegin(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 50674ae819SStefano Zampini ierr = PCBDDCScatterCoarseDataEnd(pc,local_primal_vec,pcbddc->coarse_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 51fdc635d7SStefano Zampini if (coarse_mat) { 52674ae819SStefano Zampini if (pcbddc->dbg_flag) { 53fdc635d7SStefano Zampini ierr = MatMult(coarse_mat,pcbddc->coarse_vec,pcbddc->coarse_rhs);CHKERRQ(ierr); 54674ae819SStefano Zampini ierr = VecNorm(pcbddc->coarse_rhs,NORM_2,&test_null);CHKERRQ(ierr); 55674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Vec %d coarse null space error % 1.14e\n",i,test_null);CHKERRQ(ierr); 56674ae819SStefano Zampini } 57674ae819SStefano Zampini ierr = VecCopy(pcbddc->coarse_vec,coarse_nsp_vecs[coarse_nsp_size]);CHKERRQ(ierr); 58674ae819SStefano Zampini coarse_nsp_size++; 59674ae819SStefano Zampini } 60674ae819SStefano Zampini } 61674ae819SStefano Zampini if (coarse_nsp_size > 0) { 629a7d3425SStefano Zampini ierr = PCBDDCOrthonormalizeVecs(coarse_nsp_size,coarse_nsp_vecs);CHKERRQ(ierr); 63fdc635d7SStefano Zampini ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject)coarse_mat),PETSC_FALSE,coarse_nsp_size,coarse_nsp_vecs,&tempCoarseNullSpace);CHKERRQ(ierr); 64674ae819SStefano Zampini for (i=0;i<nsp_size+1;i++) { 65674ae819SStefano Zampini ierr = VecDestroy(&coarse_nsp_vecs[i]);CHKERRQ(ierr); 66674ae819SStefano Zampini } 67674ae819SStefano Zampini } 68674ae819SStefano Zampini ierr = PetscFree(coarse_nsp_vecs);CHKERRQ(ierr); 69674ae819SStefano Zampini ierr = VecDestroy(&local_vec);CHKERRQ(ierr); 70674ae819SStefano Zampini ierr = VecDestroy(&local_primal_vec);CHKERRQ(ierr); 71674ae819SStefano Zampini *CoarseNullSpace = tempCoarseNullSpace; 72674ae819SStefano Zampini PetscFunctionReturn(0); 73674ae819SStefano Zampini } 74674ae819SStefano Zampini 75674ae819SStefano Zampini 76674ae819SStefano Zampini #undef __FUNCT__ 77674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplyNullSpaceCorrectionPC" 78674ae819SStefano Zampini static PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC pc,Vec x,Vec y) 79674ae819SStefano Zampini { 80674ae819SStefano Zampini NullSpaceCorrection_ctx pc_ctx; 81674ae819SStefano Zampini PetscErrorCode ierr; 82674ae819SStefano Zampini 83674ae819SStefano Zampini PetscFunctionBegin; 84674ae819SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 85674ae819SStefano Zampini /* E */ 86674ae819SStefano Zampini ierr = MatMultTranspose(pc_ctx->Lbasis_mat,x,pc_ctx->work_small_2);CHKERRQ(ierr); 87674ae819SStefano Zampini ierr = MatMultAdd(pc_ctx->Kbasis_mat,pc_ctx->work_small_2,x,pc_ctx->work_full_1);CHKERRQ(ierr); 88674ae819SStefano Zampini /* P^-1 */ 89674ae819SStefano Zampini ierr = PCApply(pc_ctx->local_pc,pc_ctx->work_full_1,pc_ctx->work_full_2);CHKERRQ(ierr); 90674ae819SStefano Zampini /* E^T */ 91674ae819SStefano Zampini ierr = MatMultTranspose(pc_ctx->Kbasis_mat,pc_ctx->work_full_2,pc_ctx->work_small_1);CHKERRQ(ierr); 92674ae819SStefano Zampini ierr = VecScale(pc_ctx->work_small_1,-1.0);CHKERRQ(ierr); 93674ae819SStefano Zampini ierr = MatMultAdd(pc_ctx->Lbasis_mat,pc_ctx->work_small_1,pc_ctx->work_full_2,pc_ctx->work_full_1);CHKERRQ(ierr); 94674ae819SStefano Zampini /* Sum contributions */ 95674ae819SStefano Zampini ierr = MatMultAdd(pc_ctx->basis_mat,pc_ctx->work_small_2,pc_ctx->work_full_1,y);CHKERRQ(ierr); 96674ae819SStefano Zampini PetscFunctionReturn(0); 97674ae819SStefano Zampini } 98674ae819SStefano Zampini 99674ae819SStefano Zampini #undef __FUNCT__ 100674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyNullSpaceCorrectionPC" 101674ae819SStefano Zampini static PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC pc) 102674ae819SStefano Zampini { 103674ae819SStefano Zampini NullSpaceCorrection_ctx pc_ctx; 104674ae819SStefano Zampini PetscErrorCode ierr; 105674ae819SStefano Zampini 106674ae819SStefano Zampini PetscFunctionBegin; 107674ae819SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 108674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_small_1);CHKERRQ(ierr); 109674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_small_2);CHKERRQ(ierr); 110674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_full_1);CHKERRQ(ierr); 111674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->work_full_2);CHKERRQ(ierr); 112674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->basis_mat);CHKERRQ(ierr); 113674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->Lbasis_mat);CHKERRQ(ierr); 114674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->Kbasis_mat);CHKERRQ(ierr); 115674ae819SStefano Zampini ierr = PCDestroy(&pc_ctx->local_pc);CHKERRQ(ierr); 116674ae819SStefano Zampini ierr = PetscFree(pc_ctx);CHKERRQ(ierr); 117674ae819SStefano Zampini PetscFunctionReturn(0); 118674ae819SStefano Zampini } 119674ae819SStefano Zampini 12028509bceSStefano Zampini /* 12128509bceSStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCApplyNullSpaceCorrectionPC(PC,Vec,Vec); 12228509bceSStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCDestroyNullSpaceCorrectionPC(PC); 12328509bceSStefano Zampini */ 124674ae819SStefano Zampini 125674ae819SStefano Zampini #undef __FUNCT__ 126674ae819SStefano Zampini #define __FUNCT__ "PCBDDCNullSpaceAssembleCorrection" 127674ae819SStefano Zampini PetscErrorCode PCBDDCNullSpaceAssembleCorrection(PC pc,IS local_dofs) 128674ae819SStefano Zampini { 129674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 130674ae819SStefano Zampini PC_IS *pcis = (PC_IS*)pc->data; 131674ae819SStefano Zampini Mat_IS* matis = (Mat_IS*)pc->pmat->data; 132674ae819SStefano Zampini KSP *local_ksp; 133674ae819SStefano Zampini PC newpc; 134674ae819SStefano Zampini NullSpaceCorrection_ctx shell_ctx; 135674ae819SStefano Zampini Mat local_mat,local_pmat,small_mat,inv_small_mat; 136674ae819SStefano Zampini MatStructure local_mat_struct; 137674ae819SStefano Zampini Vec work1,work2; 138674ae819SStefano Zampini const Vec *nullvecs; 139674ae819SStefano Zampini VecScatter scatter_ctx; 140674ae819SStefano Zampini IS is_aux; 141674ae819SStefano Zampini MatFactorInfo matinfo; 142674ae819SStefano Zampini PetscScalar *basis_mat,*Kbasis_mat,*array,*array_mat; 143674ae819SStefano Zampini PetscScalar one = 1.0,zero = 0.0, m_one = -1.0; 144674ae819SStefano Zampini PetscInt basis_dofs,basis_size,nnsp_size,i,k,n_I,n_R; 145674ae819SStefano Zampini PetscBool nnsp_has_cnst; 146674ae819SStefano Zampini PetscErrorCode ierr; 147674ae819SStefano Zampini 148674ae819SStefano Zampini PetscFunctionBegin; 149674ae819SStefano Zampini /* Infer the local solver */ 150674ae819SStefano Zampini ierr = ISGetSize(local_dofs,&basis_dofs);CHKERRQ(ierr); 151674ae819SStefano Zampini ierr = VecGetSize(pcis->vec1_D,&n_I);CHKERRQ(ierr); 152674ae819SStefano Zampini ierr = VecGetSize(pcbddc->vec1_R,&n_R);CHKERRQ(ierr); 153674ae819SStefano Zampini if (basis_dofs == n_I) { 154674ae819SStefano Zampini /* Dirichlet solver */ 155674ae819SStefano Zampini local_ksp = &pcbddc->ksp_D; 156674ae819SStefano Zampini } else if (basis_dofs == n_R) { 157674ae819SStefano Zampini /* Neumann solver */ 158674ae819SStefano Zampini local_ksp = &pcbddc->ksp_R; 159674ae819SStefano Zampini } else { 160674ae819SStefano 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); 161674ae819SStefano Zampini } 162674ae819SStefano Zampini ierr = KSPGetOperators(*local_ksp,&local_mat,&local_pmat,&local_mat_struct);CHKERRQ(ierr); 163674ae819SStefano Zampini 164674ae819SStefano Zampini /* Get null space vecs */ 165674ae819SStefano Zampini ierr = MatNullSpaceGetVecs(pcbddc->NullSpace,&nnsp_has_cnst,&nnsp_size,&nullvecs);CHKERRQ(ierr); 166674ae819SStefano Zampini basis_size = nnsp_size; 167674ae819SStefano Zampini if (nnsp_has_cnst) { 168674ae819SStefano Zampini basis_size++; 169674ae819SStefano Zampini } 170674ae819SStefano Zampini 171b8ffe317SStefano Zampini if (basis_dofs) { 172674ae819SStefano Zampini /* Create shell ctx */ 173674ae819SStefano Zampini ierr = PetscMalloc(sizeof(*shell_ctx),&shell_ctx);CHKERRQ(ierr); 174674ae819SStefano Zampini 175674ae819SStefano Zampini /* Create work vectors in shell context */ 176674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_small_1);CHKERRQ(ierr); 177674ae819SStefano Zampini ierr = VecSetSizes(shell_ctx->work_small_1,basis_size,basis_size);CHKERRQ(ierr); 178674ae819SStefano Zampini ierr = VecSetType(shell_ctx->work_small_1,VECSEQ);CHKERRQ(ierr); 179674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_small_1,&shell_ctx->work_small_2);CHKERRQ(ierr); 180674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&shell_ctx->work_full_1);CHKERRQ(ierr); 181674ae819SStefano Zampini ierr = VecSetSizes(shell_ctx->work_full_1,basis_dofs,basis_dofs);CHKERRQ(ierr); 182674ae819SStefano Zampini ierr = VecSetType(shell_ctx->work_full_1,VECSEQ);CHKERRQ(ierr); 183674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&shell_ctx->work_full_2);CHKERRQ(ierr); 184674ae819SStefano Zampini 185674ae819SStefano Zampini /* Allocate workspace */ 186674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->basis_mat );CHKERRQ(ierr); 187674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_dofs,basis_size,NULL,&shell_ctx->Kbasis_mat);CHKERRQ(ierr); 188674ae819SStefano Zampini ierr = MatDenseGetArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr); 189674ae819SStefano Zampini ierr = MatDenseGetArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr); 190674ae819SStefano Zampini 191674ae819SStefano Zampini /* Restrict local null space on selected dofs (Dirichlet or Neumann) 192674ae819SStefano Zampini and compute matrices N and K*N */ 193674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work1);CHKERRQ(ierr); 194674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work2);CHKERRQ(ierr); 195674ae819SStefano Zampini ierr = VecScatterCreate(pcis->vec1_N,local_dofs,work1,(IS)0,&scatter_ctx);CHKERRQ(ierr); 196b8ffe317SStefano Zampini } 197b8ffe317SStefano Zampini 198674ae819SStefano Zampini for (k=0;k<nnsp_size;k++) { 199674ae819SStefano Zampini ierr = VecScatterBegin(matis->ctx,nullvecs[k],pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 200674ae819SStefano Zampini ierr = VecScatterEnd(matis->ctx,nullvecs[k],pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 201b8ffe317SStefano Zampini if (basis_dofs) { 202674ae819SStefano Zampini ierr = VecPlaceArray(work1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr); 203674ae819SStefano Zampini ierr = VecScatterBegin(scatter_ctx,pcis->vec1_N,work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 204674ae819SStefano Zampini ierr = VecScatterEnd(scatter_ctx,pcis->vec1_N,work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 205674ae819SStefano Zampini ierr = VecPlaceArray(work2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr); 206674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr); 207674ae819SStefano Zampini ierr = VecResetArray(work1);CHKERRQ(ierr); 208674ae819SStefano Zampini ierr = VecResetArray(work2);CHKERRQ(ierr); 209674ae819SStefano Zampini } 210b8ffe317SStefano Zampini } 211b8ffe317SStefano Zampini 212b8ffe317SStefano Zampini if (basis_dofs) { 213674ae819SStefano Zampini if (nnsp_has_cnst) { 214674ae819SStefano Zampini ierr = VecPlaceArray(work1,(const PetscScalar*)&basis_mat[k*basis_dofs]);CHKERRQ(ierr); 215674ae819SStefano Zampini ierr = VecSet(work1,one);CHKERRQ(ierr); 216674ae819SStefano Zampini ierr = VecPlaceArray(work2,(const PetscScalar*)&Kbasis_mat[k*basis_dofs]);CHKERRQ(ierr); 217674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr); 218674ae819SStefano Zampini ierr = VecResetArray(work1);CHKERRQ(ierr); 219674ae819SStefano Zampini ierr = VecResetArray(work2);CHKERRQ(ierr); 220674ae819SStefano Zampini } 221674ae819SStefano Zampini ierr = VecDestroy(&work1);CHKERRQ(ierr); 222674ae819SStefano Zampini ierr = VecDestroy(&work2);CHKERRQ(ierr); 223674ae819SStefano Zampini ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr); 224674ae819SStefano Zampini ierr = MatDenseRestoreArray(shell_ctx->basis_mat,&basis_mat);CHKERRQ(ierr); 225674ae819SStefano Zampini ierr = MatDenseRestoreArray(shell_ctx->Kbasis_mat,&Kbasis_mat);CHKERRQ(ierr); 226674ae819SStefano Zampini 227674ae819SStefano Zampini /* Assemble another Mat object in shell context */ 228674ae819SStefano Zampini ierr = MatTransposeMatMult(shell_ctx->basis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&small_mat);CHKERRQ(ierr); 229674ae819SStefano Zampini ierr = MatFactorInfoInitialize(&matinfo);CHKERRQ(ierr); 230674ae819SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,basis_size,0,1,&is_aux);CHKERRQ(ierr); 231674ae819SStefano Zampini ierr = MatLUFactor(small_mat,is_aux,is_aux,&matinfo);CHKERRQ(ierr); 232674ae819SStefano Zampini ierr = ISDestroy(&is_aux);CHKERRQ(ierr); 233*785e854fSJed Brown ierr = PetscMalloc1(basis_size*basis_size,&array_mat);CHKERRQ(ierr); 234674ae819SStefano Zampini for (k=0;k<basis_size;k++) { 235674ae819SStefano Zampini ierr = VecSet(shell_ctx->work_small_1,zero);CHKERRQ(ierr); 236674ae819SStefano Zampini ierr = VecSetValue(shell_ctx->work_small_1,k,one,INSERT_VALUES);CHKERRQ(ierr); 237674ae819SStefano Zampini ierr = VecAssemblyBegin(shell_ctx->work_small_1);CHKERRQ(ierr); 238674ae819SStefano Zampini ierr = VecAssemblyEnd(shell_ctx->work_small_1);CHKERRQ(ierr); 239674ae819SStefano Zampini ierr = MatSolve(small_mat,shell_ctx->work_small_1,shell_ctx->work_small_2);CHKERRQ(ierr); 240674ae819SStefano Zampini ierr = VecGetArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr); 241674ae819SStefano Zampini for (i=0;i<basis_size;i++) { 242674ae819SStefano Zampini array_mat[i*basis_size+k]=array[i]; 243674ae819SStefano Zampini } 244674ae819SStefano Zampini ierr = VecRestoreArrayRead(shell_ctx->work_small_2,(const PetscScalar**)&array);CHKERRQ(ierr); 245674ae819SStefano Zampini } 246674ae819SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,basis_size,basis_size,array_mat,&inv_small_mat);CHKERRQ(ierr); 247674ae819SStefano Zampini ierr = MatMatMult(shell_ctx->basis_mat,inv_small_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&shell_ctx->Lbasis_mat);CHKERRQ(ierr); 248674ae819SStefano Zampini ierr = PetscFree(array_mat);CHKERRQ(ierr); 249674ae819SStefano Zampini ierr = MatDestroy(&inv_small_mat);CHKERRQ(ierr); 250674ae819SStefano Zampini ierr = MatDestroy(&small_mat);CHKERRQ(ierr); 251674ae819SStefano Zampini ierr = MatScale(shell_ctx->Kbasis_mat,m_one);CHKERRQ(ierr); 252674ae819SStefano Zampini 253674ae819SStefano Zampini /* Rebuild local PC */ 254674ae819SStefano Zampini ierr = KSPGetPC(*local_ksp,&shell_ctx->local_pc);CHKERRQ(ierr); 255674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)shell_ctx->local_pc);CHKERRQ(ierr); 256674ae819SStefano Zampini ierr = PCCreate(PETSC_COMM_SELF,&newpc);CHKERRQ(ierr); 257674ae819SStefano Zampini ierr = PCSetOperators(newpc,local_mat,local_mat,SAME_PRECONDITIONER);CHKERRQ(ierr); 258674ae819SStefano Zampini ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 259674ae819SStefano Zampini ierr = PCShellSetContext(newpc,shell_ctx);CHKERRQ(ierr); 260674ae819SStefano Zampini ierr = PCShellSetApply(newpc,PCBDDCApplyNullSpaceCorrectionPC);CHKERRQ(ierr); 261674ae819SStefano Zampini ierr = PCShellSetDestroy(newpc,PCBDDCDestroyNullSpaceCorrectionPC);CHKERRQ(ierr); 262674ae819SStefano Zampini ierr = PCSetUp(newpc);CHKERRQ(ierr); 263674ae819SStefano Zampini ierr = KSPSetPC(*local_ksp,newpc);CHKERRQ(ierr); 264674ae819SStefano Zampini ierr = PCDestroy(&newpc);CHKERRQ(ierr); 265674ae819SStefano Zampini ierr = KSPSetUp(*local_ksp);CHKERRQ(ierr); 266b8ffe317SStefano Zampini } 267674ae819SStefano Zampini /* test */ 268b8ffe317SStefano Zampini if (pcbddc->dbg_flag && basis_dofs) { 269674ae819SStefano Zampini KSP check_ksp; 270674ae819SStefano Zampini PC check_pc; 271674ae819SStefano Zampini Mat test_mat; 272674ae819SStefano Zampini Vec work3; 273674ae819SStefano Zampini PetscReal test_err,lambda_min,lambda_max; 274674ae819SStefano Zampini PetscBool setsym,issym=PETSC_FALSE; 275b8ffe317SStefano Zampini PetscInt tabs; 276674ae819SStefano Zampini 277b8ffe317SStefano Zampini ierr = PetscViewerASCIIGetTab(pcbddc->dbg_viewer,&tabs);CHKERRQ(ierr); 278674ae819SStefano Zampini ierr = KSPGetPC(*local_ksp,&check_pc);CHKERRQ(ierr); 279674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work1);CHKERRQ(ierr); 280674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work2);CHKERRQ(ierr); 281674ae819SStefano Zampini ierr = VecDuplicate(shell_ctx->work_full_1,&work3);CHKERRQ(ierr); 282674ae819SStefano Zampini ierr = VecSetRandom(shell_ctx->work_small_1,NULL);CHKERRQ(ierr); 283674ae819SStefano Zampini ierr = MatMult(shell_ctx->basis_mat,shell_ctx->work_small_1,work1);CHKERRQ(ierr); 284674ae819SStefano Zampini ierr = VecCopy(work1,work2);CHKERRQ(ierr); 285674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work3);CHKERRQ(ierr); 286674ae819SStefano Zampini ierr = PCApply(check_pc,work3,work1);CHKERRQ(ierr); 287674ae819SStefano Zampini ierr = VecAXPY(work1,m_one,work2);CHKERRQ(ierr); 288674ae819SStefano Zampini ierr = VecNorm(work1,NORM_INFINITY,&test_err);CHKERRQ(ierr); 289b8ffe317SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d error for nullspace correction for ",PetscGlobalRank); 290b8ffe317SStefano Zampini ierr = PetscViewerASCIIUseTabs(pcbddc->dbg_viewer,PETSC_FALSE);CHKERRQ(ierr); 291674ae819SStefano Zampini if (basis_dofs == n_I) { 292b8ffe317SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Dirichlet "); 293674ae819SStefano Zampini } else { 294b8ffe317SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Neumann "); 295674ae819SStefano Zampini } 296b8ffe317SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"solver is :%1.14e\n",test_err); 297b8ffe317SStefano Zampini ierr = PetscViewerASCIISetTab(pcbddc->dbg_viewer,tabs);CHKERRQ(ierr); 298b8ffe317SStefano Zampini ierr = PetscViewerASCIIUseTabs(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 299674ae819SStefano Zampini 300674ae819SStefano Zampini ierr = MatTransposeMatMult(shell_ctx->Lbasis_mat,shell_ctx->Kbasis_mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&test_mat);CHKERRQ(ierr); 301674ae819SStefano Zampini ierr = MatShift(test_mat,one);CHKERRQ(ierr); 302674ae819SStefano Zampini ierr = MatNorm(test_mat,NORM_INFINITY,&test_err);CHKERRQ(ierr); 303674ae819SStefano Zampini ierr = MatDestroy(&test_mat);CHKERRQ(ierr); 304b8ffe317SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d error for nullspace matrices is :%1.14e\n",PetscGlobalRank,test_err); 305674ae819SStefano Zampini 306674ae819SStefano Zampini /* Create ksp object suitable for extreme eigenvalues' estimation */ 307674ae819SStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&check_ksp);CHKERRQ(ierr); 308674ae819SStefano Zampini ierr = KSPSetOperators(check_ksp,local_mat,local_mat,SAME_PRECONDITIONER);CHKERRQ(ierr); 309674ae819SStefano Zampini ierr = KSPSetTolerances(check_ksp,1.e-8,1.e-8,PETSC_DEFAULT,basis_dofs);CHKERRQ(ierr); 310674ae819SStefano Zampini ierr = KSPSetComputeSingularValues(check_ksp,PETSC_TRUE);CHKERRQ(ierr); 311674ae819SStefano Zampini ierr = MatIsSymmetricKnown(pc->pmat,&setsym,&issym);CHKERRQ(ierr); 312674ae819SStefano Zampini if (issym) { 313674ae819SStefano Zampini ierr = KSPSetType(check_ksp,KSPCG);CHKERRQ(ierr); 314674ae819SStefano Zampini } 315674ae819SStefano Zampini ierr = KSPSetPC(check_ksp,check_pc);CHKERRQ(ierr); 316674ae819SStefano Zampini ierr = KSPSetUp(check_ksp);CHKERRQ(ierr); 317674ae819SStefano Zampini ierr = VecSetRandom(work1,NULL);CHKERRQ(ierr); 318674ae819SStefano Zampini ierr = MatMult(local_mat,work1,work2);CHKERRQ(ierr); 319674ae819SStefano Zampini ierr = KSPSolve(check_ksp,work2,work2);CHKERRQ(ierr); 320674ae819SStefano Zampini ierr = VecAXPY(work2,m_one,work1);CHKERRQ(ierr); 321674ae819SStefano Zampini ierr = VecNorm(work2,NORM_INFINITY,&test_err);CHKERRQ(ierr); 322674ae819SStefano Zampini ierr = KSPComputeExtremeSingularValues(check_ksp,&lambda_max,&lambda_min);CHKERRQ(ierr); 323674ae819SStefano Zampini ierr = KSPGetIterationNumber(check_ksp,&k);CHKERRQ(ierr); 324b8ffe317SStefano 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); 325674ae819SStefano Zampini ierr = KSPDestroy(&check_ksp);CHKERRQ(ierr); 326674ae819SStefano Zampini ierr = VecDestroy(&work1);CHKERRQ(ierr); 327674ae819SStefano Zampini ierr = VecDestroy(&work2);CHKERRQ(ierr); 328674ae819SStefano Zampini ierr = VecDestroy(&work3);CHKERRQ(ierr); 329674ae819SStefano Zampini } 330b8ffe317SStefano Zampini /* all processes shoud call this, even the void ones */ 331b8ffe317SStefano Zampini if (pcbddc->dbg_flag) { 332b8ffe317SStefano Zampini ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 333b8ffe317SStefano Zampini } 334674ae819SStefano Zampini PetscFunctionReturn(0); 335674ae819SStefano Zampini } 336674ae819SStefano Zampini 337674ae819SStefano Zampini #undef __FUNCT__ 338674ae819SStefano Zampini #define __FUNCT__ "PCBDDCNullSpaceAdaptGlobal" 339674ae819SStefano Zampini PetscErrorCode PCBDDCNullSpaceAdaptGlobal(PC pc) 340674ae819SStefano Zampini { 341674ae819SStefano Zampini PC_IS* pcis = (PC_IS*) (pc->data); 342674ae819SStefano Zampini PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 343674ae819SStefano Zampini KSP inv_change; 344674ae819SStefano Zampini PC pc_change; 345674ae819SStefano Zampini const Vec *nsp_vecs; 346674ae819SStefano Zampini Vec *new_nsp_vecs; 347674ae819SStefano Zampini PetscInt i,nsp_size,new_nsp_size,start_new; 348674ae819SStefano Zampini PetscBool nsp_has_cnst; 349674ae819SStefano Zampini MatNullSpace new_nsp; 350674ae819SStefano Zampini PetscErrorCode ierr; 351674ae819SStefano Zampini MPI_Comm comm; 352674ae819SStefano Zampini 353674ae819SStefano Zampini PetscFunctionBegin; 354674ae819SStefano Zampini /* create KSP for change of basis */ 355674ae819SStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&inv_change);CHKERRQ(ierr); 356674ae819SStefano Zampini ierr = KSPSetOperators(inv_change,pcbddc->ChangeOfBasisMatrix,pcbddc->ChangeOfBasisMatrix,SAME_PRECONDITIONER);CHKERRQ(ierr); 357674ae819SStefano Zampini ierr = KSPSetType(inv_change,KSPPREONLY);CHKERRQ(ierr); 358674ae819SStefano Zampini ierr = KSPGetPC(inv_change,&pc_change);CHKERRQ(ierr); 359674ae819SStefano Zampini ierr = PCSetType(pc_change,PCLU);CHKERRQ(ierr); 360674ae819SStefano Zampini ierr = KSPSetUp(inv_change);CHKERRQ(ierr); 361674ae819SStefano Zampini /* get nullspace and transform it */ 362674ae819SStefano Zampini ierr = MatNullSpaceGetVecs(pcbddc->NullSpace,&nsp_has_cnst,&nsp_size,&nsp_vecs);CHKERRQ(ierr); 363674ae819SStefano Zampini new_nsp_size = nsp_size; 364674ae819SStefano Zampini if (nsp_has_cnst) { 365674ae819SStefano Zampini new_nsp_size++; 366674ae819SStefano Zampini } 367*785e854fSJed Brown ierr = PetscMalloc1(new_nsp_size,&new_nsp_vecs);CHKERRQ(ierr); 368674ae819SStefano Zampini for (i=0;i<new_nsp_size;i++) { 369674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&new_nsp_vecs[i]);CHKERRQ(ierr); 370674ae819SStefano Zampini } 371674ae819SStefano Zampini start_new = 0; 372674ae819SStefano Zampini if (nsp_has_cnst) { 373674ae819SStefano Zampini start_new = 1; 374674ae819SStefano Zampini ierr = VecSet(new_nsp_vecs[0],1.0);CHKERRQ(ierr); 375674ae819SStefano Zampini ierr = VecSet(pcis->vec1_B,1.0);CHKERRQ(ierr); 376674ae819SStefano Zampini ierr = KSPSolve(inv_change,pcis->vec1_B,pcis->vec1_B); 377674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[0],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 378674ae819SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[0],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 379674ae819SStefano Zampini } 380674ae819SStefano Zampini for (i=0;i<nsp_size;i++) { 381674ae819SStefano Zampini ierr = VecCopy(nsp_vecs[i],new_nsp_vecs[i+start_new]);CHKERRQ(ierr); 382674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,nsp_vecs[i],pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 383674ae819SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,nsp_vecs[i],pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 384674ae819SStefano Zampini ierr = KSPSolve(inv_change,pcis->vec1_B,pcis->vec1_B); 385674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[i+start_new],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 386674ae819SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,new_nsp_vecs[i+start_new],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 387674ae819SStefano Zampini } 3889a7d3425SStefano Zampini ierr = PCBDDCOrthonormalizeVecs(new_nsp_size,new_nsp_vecs);CHKERRQ(ierr); 389674ae819SStefano Zampini #if 0 390674ae819SStefano Zampini PetscBool nsp_t=PETSC_FALSE; 391674ae819SStefano Zampini ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr); 392674ae819SStefano Zampini printf("Original Null Space test: %d\n",nsp_t); 393674ae819SStefano Zampini Mat temp_mat; 394674ae819SStefano Zampini Mat_IS* matis = (Mat_IS*)pc->pmat->data; 395674ae819SStefano Zampini temp_mat = matis->A; 396674ae819SStefano Zampini matis->A = pcbddc->local_mat; 397674ae819SStefano Zampini pcbddc->local_mat = temp_mat; 398674ae819SStefano Zampini ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr); 399674ae819SStefano Zampini printf("Original Null Space, mat changed test: %d\n",nsp_t); 400674ae819SStefano Zampini { 401674ae819SStefano Zampini PetscReal test_norm; 402674ae819SStefano Zampini for (i=0;i<new_nsp_size;i++) { 403674ae819SStefano Zampini ierr = MatMult(pc->pmat,new_nsp_vecs[i],pcis->vec1_global);CHKERRQ(ierr); 404674ae819SStefano Zampini ierr = VecNorm(pcis->vec1_global,NORM_2,&test_norm);CHKERRQ(ierr); 405674ae819SStefano Zampini if (test_norm > 1.e-12) { 406674ae819SStefano Zampini printf("------------ERROR VEC %d------------------\n",i); 407674ae819SStefano Zampini ierr = VecView(pcis->vec1_global,PETSC_VIEWER_STDOUT_WORLD); 408674ae819SStefano Zampini printf("------------------------------------------\n"); 409674ae819SStefano Zampini } 410674ae819SStefano Zampini } 411674ae819SStefano Zampini } 412674ae819SStefano Zampini #endif 413674ae819SStefano Zampini 414674ae819SStefano Zampini ierr = KSPDestroy(&inv_change);CHKERRQ(ierr); 415674ae819SStefano Zampini ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 416674ae819SStefano Zampini ierr = MatNullSpaceCreate(comm,PETSC_FALSE,new_nsp_size,new_nsp_vecs,&new_nsp);CHKERRQ(ierr); 417674ae819SStefano Zampini ierr = PCBDDCSetNullSpace(pc,new_nsp);CHKERRQ(ierr); 418674ae819SStefano Zampini ierr = MatNullSpaceDestroy(&new_nsp);CHKERRQ(ierr); 419674ae819SStefano Zampini #if 0 420674ae819SStefano Zampini ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr); 421674ae819SStefano Zampini printf("New Null Space, mat changed: %d\n",nsp_t); 422674ae819SStefano Zampini temp_mat = matis->A; 423674ae819SStefano Zampini matis->A = pcbddc->local_mat; 424674ae819SStefano Zampini pcbddc->local_mat = temp_mat; 425674ae819SStefano Zampini ierr = MatNullSpaceTest(pcbddc->NullSpace,pc->pmat,&nsp_t);CHKERRQ(ierr); 426674ae819SStefano Zampini printf("New Null Space, mat original: %d\n",nsp_t); 427674ae819SStefano Zampini #endif 428674ae819SStefano Zampini 429674ae819SStefano Zampini for (i=0;i<new_nsp_size;i++) { 430674ae819SStefano Zampini ierr = VecDestroy(&new_nsp_vecs[i]);CHKERRQ(ierr); 431674ae819SStefano Zampini } 432674ae819SStefano Zampini ierr = PetscFree(new_nsp_vecs);CHKERRQ(ierr); 433674ae819SStefano Zampini PetscFunctionReturn(0); 434674ae819SStefano Zampini } 435