1 #include "bddc.h" 2 #include "bddcprivate.h" 3 #include <petscblaslapack.h> 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "PCBDDCResetCustomization" 7 PetscErrorCode PCBDDCResetCustomization(PC pc) 8 { 9 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 10 PetscInt i; 11 PetscErrorCode ierr; 12 13 PetscFunctionBegin; 14 ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr); 15 ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 16 ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 17 ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 18 ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 19 for (i=0;i<pcbddc->n_ISForDofs;i++) { 20 ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 21 } 22 ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 23 PetscFunctionReturn(0); 24 } 25 26 #undef __FUNCT__ 27 #define __FUNCT__ "PCBDDCResetTopography" 28 PetscErrorCode PCBDDCResetTopography(PC pc) 29 { 30 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 31 PetscErrorCode ierr; 32 33 PetscFunctionBegin; 34 ierr = MatDestroy(&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr); 35 ierr = MatDestroy(&pcbddc->ConstraintMatrix);CHKERRQ(ierr); 36 ierr = PCBDDCGraphReset(pcbddc->mat_graph);CHKERRQ(ierr); 37 PetscFunctionReturn(0); 38 } 39 40 #undef __FUNCT__ 41 #define __FUNCT__ "PCBDDCResetSolvers" 42 PetscErrorCode PCBDDCResetSolvers(PC pc) 43 { 44 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 45 PetscErrorCode ierr; 46 47 PetscFunctionBegin; 48 ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 49 ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 50 ierr = VecDestroy(&pcbddc->coarse_vec);CHKERRQ(ierr); 51 ierr = VecDestroy(&pcbddc->coarse_rhs);CHKERRQ(ierr); 52 ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr); 53 ierr = MatDestroy(&pcbddc->coarse_mat);CHKERRQ(ierr); 54 ierr = MatDestroy(&pcbddc->coarse_phi_B);CHKERRQ(ierr); 55 ierr = MatDestroy(&pcbddc->coarse_phi_D);CHKERRQ(ierr); 56 ierr = MatDestroy(&pcbddc->coarse_psi_B);CHKERRQ(ierr); 57 ierr = MatDestroy(&pcbddc->coarse_psi_D);CHKERRQ(ierr); 58 ierr = VecDestroy(&pcbddc->vec1_P);CHKERRQ(ierr); 59 ierr = VecDestroy(&pcbddc->vec1_C);CHKERRQ(ierr); 60 ierr = MatDestroy(&pcbddc->local_auxmat1);CHKERRQ(ierr); 61 ierr = MatDestroy(&pcbddc->local_auxmat2);CHKERRQ(ierr); 62 ierr = VecDestroy(&pcbddc->vec1_R);CHKERRQ(ierr); 63 ierr = VecDestroy(&pcbddc->vec2_R);CHKERRQ(ierr); 64 ierr = VecDestroy(&pcbddc->vec4_D);CHKERRQ(ierr); 65 ierr = VecScatterDestroy(&pcbddc->R_to_B);CHKERRQ(ierr); 66 ierr = VecScatterDestroy(&pcbddc->R_to_D);CHKERRQ(ierr); 67 ierr = VecScatterDestroy(&pcbddc->coarse_loc_to_glob);CHKERRQ(ierr); 68 ierr = PetscFree(pcbddc->local_primal_indices);CHKERRQ(ierr); 69 ierr = PetscFree(pcbddc->replicated_local_primal_indices);CHKERRQ(ierr); 70 ierr = PetscFree(pcbddc->replicated_local_primal_values);CHKERRQ(ierr); 71 ierr = PetscFree(pcbddc->local_primal_displacements);CHKERRQ(ierr); 72 ierr = PetscFree(pcbddc->local_primal_sizes);CHKERRQ(ierr); 73 PetscFunctionReturn(0); 74 } 75 76 #undef __FUNCT__ 77 #define __FUNCT__ "PCBDDCCreateWorkVectors" 78 PetscErrorCode PCBDDCCreateWorkVectors(PC pc) 79 { 80 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 81 PC_IS *pcis = (PC_IS*)pc->data; 82 VecType impVecType; 83 PetscInt n_vertices,n_constraints,local_primal_size,n_R; 84 PetscErrorCode ierr; 85 86 PetscFunctionBegin; 87 ierr = PCBDDCGetPrimalVerticesLocalIdx(pc,&n_vertices,NULL);CHKERRQ(ierr); 88 ierr = PCBDDCGetPrimalConstraintsLocalIdx(pc,&n_constraints,NULL);CHKERRQ(ierr); 89 local_primal_size = n_constraints+n_vertices; 90 n_R = pcis->n-n_vertices; 91 /* parallel work vectors used in presolve. TODO: move outside */ 92 ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 93 ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 94 /* local work vectors */ 95 ierr = VecGetType(pcis->vec1_N,&impVecType);CHKERRQ(ierr); 96 ierr = VecDuplicate(pcis->vec1_D,&pcbddc->vec4_D);CHKERRQ(ierr); 97 ierr = VecCreate(PetscObjectComm((PetscObject)pcis->vec1_N),&pcbddc->vec1_R);CHKERRQ(ierr); 98 ierr = VecSetSizes(pcbddc->vec1_R,PETSC_DECIDE,n_R);CHKERRQ(ierr); 99 ierr = VecSetType(pcbddc->vec1_R,impVecType);CHKERRQ(ierr); 100 ierr = VecDuplicate(pcbddc->vec1_R,&pcbddc->vec2_R);CHKERRQ(ierr); 101 ierr = VecCreate(PetscObjectComm((PetscObject)pcis->vec1_N),&pcbddc->vec1_P);CHKERRQ(ierr); 102 ierr = VecSetSizes(pcbddc->vec1_P,PETSC_DECIDE,local_primal_size);CHKERRQ(ierr); 103 ierr = VecSetType(pcbddc->vec1_P,impVecType);CHKERRQ(ierr); 104 if (n_constraints) { 105 ierr = VecCreate(PetscObjectComm((PetscObject)pcis->vec1_N),&pcbddc->vec1_C);CHKERRQ(ierr); 106 ierr = VecSetSizes(pcbddc->vec1_C,PETSC_DECIDE,n_constraints);CHKERRQ(ierr); 107 ierr = VecSetType(pcbddc->vec1_C,impVecType);CHKERRQ(ierr); 108 } 109 PetscFunctionReturn(0); 110 } 111 112 #undef __FUNCT__ 113 #define __FUNCT__ "PCBDDCSetUpLocalMatrices" 114 PetscErrorCode PCBDDCSetUpLocalMatrices(PC pc) 115 { 116 PC_IS* pcis = (PC_IS*)(pc->data); 117 PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 118 Mat_IS* matis = (Mat_IS*)pc->pmat->data; 119 /* manage repeated solves */ 120 MatReuse reuse; 121 MatStructure matstruct; 122 PetscErrorCode ierr; 123 124 PetscFunctionBegin; 125 /* get mat flags */ 126 ierr = PCGetOperators(pc,NULL,NULL,&matstruct);CHKERRQ(ierr); 127 reuse = MAT_INITIAL_MATRIX; 128 if (pc->setupcalled) { 129 /* when matstruct is SAME_PRECONDITIONER, we shouldn't be here */ 130 if (matstruct == SAME_PRECONDITIONER) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"This should not happen"); 131 if (matstruct == SAME_NONZERO_PATTERN) { 132 reuse = MAT_REUSE_MATRIX; 133 } else { 134 reuse = MAT_INITIAL_MATRIX; 135 } 136 } 137 if (reuse == MAT_INITIAL_MATRIX) { 138 ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr); 139 ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr); 140 ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr); 141 ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr); 142 ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 143 } 144 145 /* transform local matrices if needed */ 146 if (pcbddc->use_change_of_basis) { 147 Mat change_mat_all; 148 PetscScalar *row_cmat_values; 149 PetscInt *row_cmat_indices; 150 PetscInt *nnz,*is_indices,*temp_indices; 151 PetscInt i,j,k,n_D,n_B; 152 153 /* Get Non-overlapping dimensions */ 154 n_B = pcis->n_B; 155 n_D = pcis->n-n_B; 156 157 /* compute nonzero structure of change of basis on all local nodes */ 158 ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 159 ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 160 for (i=0;i<n_D;i++) nnz[is_indices[i]] = 1; 161 ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 162 ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 163 k=1; 164 for (i=0;i<n_B;i++) { 165 ierr = MatGetRow(pcbddc->ChangeOfBasisMatrix,i,&j,NULL,NULL);CHKERRQ(ierr); 166 nnz[is_indices[i]]=j; 167 if (k < j) k = j; 168 ierr = MatRestoreRow(pcbddc->ChangeOfBasisMatrix,i,&j,NULL,NULL);CHKERRQ(ierr); 169 } 170 ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 171 /* assemble change of basis matrix on the whole set of local dofs */ 172 ierr = PetscMalloc(k*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr); 173 ierr = MatCreate(PETSC_COMM_SELF,&change_mat_all);CHKERRQ(ierr); 174 ierr = MatSetSizes(change_mat_all,pcis->n,pcis->n,pcis->n,pcis->n);CHKERRQ(ierr); 175 ierr = MatSetType(change_mat_all,MATSEQAIJ);CHKERRQ(ierr); 176 ierr = MatSeqAIJSetPreallocation(change_mat_all,0,nnz);CHKERRQ(ierr); 177 ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 178 for (i=0;i<n_D;i++) { 179 ierr = MatSetValue(change_mat_all,is_indices[i],is_indices[i],1.0,INSERT_VALUES);CHKERRQ(ierr); 180 } 181 ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 182 ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 183 for (i=0;i<n_B;i++) { 184 ierr = MatGetRow(pcbddc->ChangeOfBasisMatrix,i,&j,(const PetscInt**)&row_cmat_indices,(const PetscScalar**)&row_cmat_values);CHKERRQ(ierr); 185 for (k=0; k<j; k++) temp_indices[k]=is_indices[row_cmat_indices[k]]; 186 ierr = MatSetValues(change_mat_all,1,&is_indices[i],j,temp_indices,row_cmat_values,INSERT_VALUES);CHKERRQ(ierr); 187 ierr = MatRestoreRow(pcbddc->ChangeOfBasisMatrix,i,&j,(const PetscInt**)&row_cmat_indices,(const PetscScalar**)&row_cmat_values);CHKERRQ(ierr); 188 } 189 ierr = MatAssemblyBegin(change_mat_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 190 ierr = MatAssemblyEnd(change_mat_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 191 /* TODO: HOW TO WORK WITH BAIJ? PtAP not provided */ 192 ierr = MatGetBlockSize(matis->A,&i);CHKERRQ(ierr); 193 if (i==1) { 194 ierr = MatPtAP(matis->A,change_mat_all,reuse,2.0,&pcbddc->local_mat);CHKERRQ(ierr); 195 } else { 196 Mat work_mat; 197 ierr = MatConvert(matis->A,MATSEQAIJ,MAT_INITIAL_MATRIX,&work_mat);CHKERRQ(ierr); 198 ierr = MatPtAP(work_mat,change_mat_all,reuse,2.0,&pcbddc->local_mat);CHKERRQ(ierr); 199 ierr = MatDestroy(&work_mat);CHKERRQ(ierr); 200 } 201 ierr = MatDestroy(&change_mat_all);CHKERRQ(ierr); 202 ierr = PetscFree(nnz);CHKERRQ(ierr); 203 ierr = PetscFree(temp_indices);CHKERRQ(ierr); 204 } else { 205 /* without change of basis, the local matrix is unchanged */ 206 if (!pcbddc->local_mat) { 207 ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 208 pcbddc->local_mat = matis->A; 209 } 210 } 211 212 /* get submatrices */ 213 ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->A_II);CHKERRQ(ierr); 214 ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr); 215 ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr); 216 ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,reuse,&pcis->A_BB);CHKERRQ(ierr); 217 PetscFunctionReturn(0); 218 } 219 220 #undef __FUNCT__ 221 #define __FUNCT__ "PCBDDCSetUpLocalScatters" 222 PetscErrorCode PCBDDCSetUpLocalScatters(PC pc,IS* is_R_local_n) 223 { 224 PC_IS* pcis = (PC_IS*)(pc->data); 225 PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 226 IS is_R_local,is_aux1,is_aux2; 227 PetscInt *vertices,*aux_array1,*aux_array2,*is_indices,*idx_R_local; 228 PetscInt n_vertices,n_constraints,i,j,n_R,n_D,n_B; 229 PetscBool *array_bool; 230 PetscErrorCode ierr; 231 232 PetscFunctionBegin; 233 /* Set Non-overlapping dimensions */ 234 n_B = pcis->n_B; n_D = pcis->n - n_B; 235 /* get vertex indices from constraint matrix */ 236 ierr = PCBDDCGetPrimalVerticesLocalIdx(pc,&n_vertices,&vertices);CHKERRQ(ierr); 237 /* Set number of constraints */ 238 n_constraints = pcbddc->local_primal_size-n_vertices; 239 /* Dohrmann's notation: dofs splitted in R (Remaining: all dofs but the vertices) and V (Vertices) */ 240 ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&array_bool);CHKERRQ(ierr); 241 for (i=0;i<pcis->n;i++) array_bool[i] = PETSC_TRUE; 242 for (i=0;i<n_vertices;i++) array_bool[vertices[i]] = PETSC_FALSE; 243 ierr = PetscMalloc((pcis->n-n_vertices)*sizeof(PetscInt),&idx_R_local);CHKERRQ(ierr); 244 for (i=0, n_R=0; i<pcis->n; i++) { 245 if (array_bool[i]) { 246 idx_R_local[n_R] = i; 247 n_R++; 248 } 249 } 250 ierr = PetscFree(vertices);CHKERRQ(ierr); 251 ierr = ISCreateGeneral(PETSC_COMM_SELF,n_R,idx_R_local,PETSC_OWN_POINTER,&is_R_local);CHKERRQ(ierr); 252 253 /* print some info if requested */ 254 if (pcbddc->dbg_flag) { 255 ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"--------------------------------------------------\n");CHKERRQ(ierr); 256 ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 257 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d local dimensions\n",PetscGlobalRank);CHKERRQ(ierr); 258 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"local_size = %d, dirichlet_size = %d, boundary_size = %d\n",pcis->n,n_D,n_B);CHKERRQ(ierr); 259 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"r_size = %d, v_size = %d, constraints = %d, local_primal_size = %d\n",n_R,n_vertices,n_constraints,pcbddc->local_primal_size);CHKERRQ(ierr); 260 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"pcbddc->n_vertices = %d, pcbddc->n_constraints = %d\n",pcbddc->n_vertices,pcbddc->n_constraints);CHKERRQ(ierr); 261 ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 262 } 263 264 /* VecScatters pcbddc->R_to_B and (optionally) pcbddc->R_to_D */ 265 ierr = PetscMalloc((pcis->n_B-n_vertices)*sizeof(PetscInt),&aux_array1);CHKERRQ(ierr); 266 ierr = PetscMalloc((pcis->n_B-n_vertices)*sizeof(PetscInt),&aux_array2);CHKERRQ(ierr); 267 ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 268 for (i=0; i<n_D; i++) array_bool[is_indices[i]] = PETSC_FALSE; 269 ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 270 for (i=0, j=0; i<n_R; i++) { 271 if (array_bool[idx_R_local[i]]) { 272 aux_array1[j] = i; 273 j++; 274 } 275 } 276 ierr = ISCreateGeneral(PETSC_COMM_SELF,j,aux_array1,PETSC_OWN_POINTER,&is_aux1);CHKERRQ(ierr); 277 ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 278 for (i=0, j=0; i<n_B; i++) { 279 if (array_bool[is_indices[i]]) { 280 aux_array2[j] = i; j++; 281 } 282 } 283 ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 284 ierr = ISCreateGeneral(PETSC_COMM_SELF,j,aux_array2,PETSC_OWN_POINTER,&is_aux2);CHKERRQ(ierr); 285 ierr = VecScatterCreate(pcbddc->vec1_R,is_aux1,pcis->vec1_B,is_aux2,&pcbddc->R_to_B);CHKERRQ(ierr); 286 ierr = ISDestroy(&is_aux1);CHKERRQ(ierr); 287 ierr = ISDestroy(&is_aux2);CHKERRQ(ierr); 288 289 if (pcbddc->inexact_prec_type || pcbddc->dbg_flag ) { 290 ierr = PetscMalloc(n_D*sizeof(PetscInt),&aux_array1);CHKERRQ(ierr); 291 for (i=0, j=0; i<n_R; i++) { 292 if (!array_bool[idx_R_local[i]]) { 293 aux_array1[j] = i; 294 j++; 295 } 296 } 297 ierr = ISCreateGeneral(PETSC_COMM_SELF,j,aux_array1,PETSC_OWN_POINTER,&is_aux1);CHKERRQ(ierr); 298 ierr = VecScatterCreate(pcbddc->vec1_R,is_aux1,pcis->vec1_D,(IS)0,&pcbddc->R_to_D);CHKERRQ(ierr); 299 ierr = ISDestroy(&is_aux1);CHKERRQ(ierr); 300 } 301 ierr = PetscFree(array_bool);CHKERRQ(ierr); 302 *is_R_local_n = is_R_local; 303 PetscFunctionReturn(0); 304 } 305 306 #undef __FUNCT__ 307 #define __FUNCT__ "PCBDDCSetUseExactDirichlet" 308 PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool use) 309 { 310 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 311 312 PetscFunctionBegin; 313 pcbddc->use_exact_dirichlet=use; 314 PetscFunctionReturn(0); 315 } 316 317 #undef __FUNCT__ 318 #define __FUNCT__ "PCBDDCSetUpLocalSolvers" 319 PetscErrorCode PCBDDCSetUpLocalSolvers(PC pc, IS is_I_local, IS is_R_local) 320 { 321 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 322 PC_IS *pcis = (PC_IS*)pc->data; 323 PC pc_temp; 324 Mat A_RR; 325 Vec vec1,vec2,vec3; 326 MatStructure matstruct; 327 PetscScalar m_one = -1.0; 328 PetscReal value; 329 PetscInt n_D,n_R,use_exact,use_exact_reduced; 330 PetscErrorCode ierr; 331 332 PetscFunctionBegin; 333 /* Creating PC contexts for local Dirichlet and Neumann problems */ 334 ierr = PCGetOperators(pc,NULL,NULL,&matstruct);CHKERRQ(ierr); 335 336 /* DIRICHLET PROBLEM */ 337 /* Matrix for Dirichlet problem is pcis->A_II */ 338 ierr = ISGetSize(is_I_local,&n_D);CHKERRQ(ierr); 339 if (!pcbddc->ksp_D) { /* create object if not yet build */ 340 ierr = KSPCreate(PETSC_COMM_SELF,&pcbddc->ksp_D);CHKERRQ(ierr); 341 ierr = PetscObjectIncrementTabLevel((PetscObject)pcbddc->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr); 342 /* default */ 343 ierr = KSPSetType(pcbddc->ksp_D,KSPPREONLY);CHKERRQ(ierr); 344 ierr = KSPSetOptionsPrefix(pcbddc->ksp_D,"dirichlet_");CHKERRQ(ierr); 345 ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr); 346 ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr); 347 ierr = PCFactorSetReuseFill(pc_temp,PETSC_TRUE);CHKERRQ(ierr); 348 } 349 ierr = KSPSetOperators(pcbddc->ksp_D,pcis->A_II,pcis->A_II,matstruct);CHKERRQ(ierr); 350 /* Allow user's customization */ 351 ierr = KSPSetFromOptions(pcbddc->ksp_D);CHKERRQ(ierr); 352 /* umfpack interface has a bug when matrix dimension is zero. TODO solve from umfpack interface */ 353 if (!n_D) { 354 ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr); 355 ierr = PCSetType(pc_temp,PCNONE);CHKERRQ(ierr); 356 } 357 /* Set Up KSP for Dirichlet problem of BDDC */ 358 ierr = KSPSetUp(pcbddc->ksp_D);CHKERRQ(ierr); 359 /* set ksp_D into pcis data */ 360 ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr); 361 ierr = PetscObjectReference((PetscObject)pcbddc->ksp_D);CHKERRQ(ierr); 362 pcis->ksp_D = pcbddc->ksp_D; 363 364 /* NEUMANN PROBLEM */ 365 /* Matrix for Neumann problem is A_RR -> we need to create it */ 366 ierr = ISGetSize(is_R_local,&n_R);CHKERRQ(ierr); 367 ierr = MatGetSubMatrix(pcbddc->local_mat,is_R_local,is_R_local,MAT_INITIAL_MATRIX,&A_RR);CHKERRQ(ierr); 368 if (!pcbddc->ksp_R) { /* create object if not yet build */ 369 ierr = KSPCreate(PETSC_COMM_SELF,&pcbddc->ksp_R);CHKERRQ(ierr); 370 ierr = PetscObjectIncrementTabLevel((PetscObject)pcbddc->ksp_R,(PetscObject)pc,1);CHKERRQ(ierr); 371 /* default */ 372 ierr = KSPSetType(pcbddc->ksp_R,KSPPREONLY);CHKERRQ(ierr); 373 ierr = KSPSetOptionsPrefix(pcbddc->ksp_R,"neumann_");CHKERRQ(ierr); 374 ierr = KSPGetPC(pcbddc->ksp_R,&pc_temp);CHKERRQ(ierr); 375 ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr); 376 ierr = PCFactorSetReuseFill(pc_temp,PETSC_TRUE);CHKERRQ(ierr); 377 } 378 ierr = KSPSetOperators(pcbddc->ksp_R,A_RR,A_RR,matstruct);CHKERRQ(ierr); 379 /* Allow user's customization */ 380 ierr = KSPSetFromOptions(pcbddc->ksp_R);CHKERRQ(ierr); 381 /* umfpack interface has a bug when matrix dimension is zero. TODO solve from umfpack interface */ 382 if (!n_R) { 383 ierr = KSPGetPC(pcbddc->ksp_R,&pc_temp);CHKERRQ(ierr); 384 ierr = PCSetType(pc_temp,PCNONE);CHKERRQ(ierr); 385 } 386 /* Set Up KSP for Neumann problem of BDDC */ 387 ierr = KSPSetUp(pcbddc->ksp_R);CHKERRQ(ierr); 388 389 /* check Dirichlet and Neumann solvers and adapt them if a nullspace correction is needed */ 390 391 /* Dirichlet */ 392 ierr = MatGetVecs(pcis->A_II,&vec1,&vec2);CHKERRQ(ierr); 393 ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr); 394 ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr); 395 ierr = MatMult(pcis->A_II,vec1,vec2);CHKERRQ(ierr); 396 ierr = KSPSolve(pcbddc->ksp_D,vec2,vec3);CHKERRQ(ierr); 397 ierr = VecAXPY(vec3,m_one,vec1);CHKERRQ(ierr); 398 ierr = VecNorm(vec3,NORM_INFINITY,&value);CHKERRQ(ierr); 399 ierr = VecDestroy(&vec1);CHKERRQ(ierr); 400 ierr = VecDestroy(&vec2);CHKERRQ(ierr); 401 ierr = VecDestroy(&vec3);CHKERRQ(ierr); 402 /* need to be adapted? */ 403 use_exact = (PetscAbsReal(value) > 1.e-4 ? 0 : 1); 404 ierr = MPI_Allreduce(&use_exact,&use_exact_reduced,1,MPIU_INT,MPI_LAND,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 405 ierr = PCBDDCSetUseExactDirichlet(pc,(PetscBool)use_exact_reduced);CHKERRQ(ierr); 406 /* print info */ 407 if (pcbddc->dbg_flag) { 408 ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 409 ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"--------------------------------------------------\n");CHKERRQ(ierr); 410 ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Checking solution of Dirichlet and Neumann problems\n");CHKERRQ(ierr); 411 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet solve = % 1.14e \n",PetscGlobalRank,value);CHKERRQ(ierr); 412 ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 413 } 414 if (n_D && pcbddc->NullSpace && !use_exact_reduced && !pcbddc->inexact_prec_type) { 415 ierr = PCBDDCNullSpaceAssembleCorrection(pc,is_I_local);CHKERRQ(ierr); 416 } 417 418 /* Neumann */ 419 ierr = MatGetVecs(A_RR,&vec1,&vec2);CHKERRQ(ierr); 420 ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr); 421 ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr); 422 ierr = MatMult(A_RR,vec1,vec2);CHKERRQ(ierr); 423 ierr = KSPSolve(pcbddc->ksp_R,vec2,vec3);CHKERRQ(ierr); 424 ierr = VecAXPY(vec3,m_one,vec1);CHKERRQ(ierr); 425 ierr = VecNorm(vec3,NORM_INFINITY,&value);CHKERRQ(ierr); 426 ierr = VecDestroy(&vec1);CHKERRQ(ierr); 427 ierr = VecDestroy(&vec2);CHKERRQ(ierr); 428 ierr = VecDestroy(&vec3);CHKERRQ(ierr); 429 /* need to be adapted? */ 430 use_exact = (PetscAbsReal(value) > 1.e-4 ? 0 : 1); 431 if (PetscAbsReal(value) > 1.e-4) use_exact = 0; 432 ierr = MPI_Allreduce(&use_exact,&use_exact_reduced,1,MPIU_INT,MPI_LAND,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 433 /* print info */ 434 if (pcbddc->dbg_flag) { 435 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Neumann solve = % 1.14e \n",PetscGlobalRank,value);CHKERRQ(ierr); 436 ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 437 } 438 if (n_R && pcbddc->NullSpace && !use_exact_reduced) { /* is it the right logic? */ 439 ierr = PCBDDCNullSpaceAssembleCorrection(pc,is_R_local);CHKERRQ(ierr); 440 } 441 442 /* free Neumann problem's matrix */ 443 ierr = MatDestroy(&A_RR);CHKERRQ(ierr); 444 PetscFunctionReturn(0); 445 } 446 447 #undef __FUNCT__ 448 #define __FUNCT__ "PCBDDCSolveSaddlePoint" 449 static PetscErrorCode PCBDDCSolveSaddlePoint(PC pc) 450 { 451 PetscErrorCode ierr; 452 PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 453 454 PetscFunctionBegin; 455 ierr = KSPSolve(pcbddc->ksp_R,pcbddc->vec1_R,pcbddc->vec2_R);CHKERRQ(ierr); 456 if (pcbddc->local_auxmat1) { 457 ierr = MatMult(pcbddc->local_auxmat1,pcbddc->vec2_R,pcbddc->vec1_C);CHKERRQ(ierr); 458 ierr = MatMultAdd(pcbddc->local_auxmat2,pcbddc->vec1_C,pcbddc->vec2_R,pcbddc->vec2_R);CHKERRQ(ierr); 459 } 460 PetscFunctionReturn(0); 461 } 462 463 #undef __FUNCT__ 464 #define __FUNCT__ "PCBDDCApplyInterfacePreconditioner" 465 PetscErrorCode PCBDDCApplyInterfacePreconditioner(PC pc) 466 { 467 PetscErrorCode ierr; 468 PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 469 PC_IS* pcis = (PC_IS*) (pc->data); 470 const PetscScalar zero = 0.0; 471 472 PetscFunctionBegin; 473 /* Application of PHI^T (or PSI^T) */ 474 if (pcbddc->coarse_psi_B) { 475 ierr = MatMultTranspose(pcbddc->coarse_psi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr); 476 if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_psi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); } 477 } else { 478 ierr = MatMultTranspose(pcbddc->coarse_phi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr); 479 if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_phi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); } 480 } 481 /* Scatter data of coarse_rhs */ 482 if (pcbddc->coarse_rhs) { ierr = VecSet(pcbddc->coarse_rhs,zero);CHKERRQ(ierr); } 483 ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 484 485 /* Local solution on R nodes */ 486 ierr = VecSet(pcbddc->vec1_R,zero);CHKERRQ(ierr); 487 ierr = VecScatterBegin(pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 488 ierr = VecScatterEnd (pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 489 if (pcbddc->inexact_prec_type) { 490 ierr = VecScatterBegin(pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 491 ierr = VecScatterEnd (pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 492 } 493 ierr = PCBDDCSolveSaddlePoint(pc);CHKERRQ(ierr); 494 ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr); 495 ierr = VecScatterBegin(pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 496 ierr = VecScatterEnd (pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 497 if (pcbddc->inexact_prec_type) { 498 ierr = VecScatterBegin(pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 499 ierr = VecScatterEnd (pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 500 } 501 502 /* Coarse solution */ 503 ierr = PCBDDCScatterCoarseDataEnd(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 504 if (pcbddc->coarse_rhs) { /* TODO remove null space when doing multilevel */ 505 ierr = KSPSolve(pcbddc->coarse_ksp,pcbddc->coarse_rhs,pcbddc->coarse_vec);CHKERRQ(ierr); 506 } 507 ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 508 ierr = PCBDDCScatterCoarseDataEnd (pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 509 510 /* Sum contributions from two levels */ 511 ierr = MatMultAdd(pcbddc->coarse_phi_B,pcbddc->vec1_P,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr); 512 if (pcbddc->inexact_prec_type) { ierr = MatMultAdd(pcbddc->coarse_phi_D,pcbddc->vec1_P,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 513 PetscFunctionReturn(0); 514 } 515 516 #undef __FUNCT__ 517 #define __FUNCT__ "PCBDDCScatterCoarseDataBegin" 518 PetscErrorCode PCBDDCScatterCoarseDataBegin(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode) 519 { 520 PetscErrorCode ierr; 521 PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 522 523 PetscFunctionBegin; 524 switch (pcbddc->coarse_communications_type) { 525 case SCATTERS_BDDC: 526 ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr); 527 break; 528 case GATHERS_BDDC: 529 break; 530 } 531 PetscFunctionReturn(0); 532 } 533 534 #undef __FUNCT__ 535 #define __FUNCT__ "PCBDDCScatterCoarseDataEnd" 536 PetscErrorCode PCBDDCScatterCoarseDataEnd(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode) 537 { 538 PetscErrorCode ierr; 539 PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 540 PetscScalar* array_to; 541 PetscScalar* array_from; 542 MPI_Comm comm; 543 PetscInt i; 544 545 PetscFunctionBegin; 546 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 547 switch (pcbddc->coarse_communications_type) { 548 case SCATTERS_BDDC: 549 ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr); 550 break; 551 case GATHERS_BDDC: 552 if (vec_from) { 553 ierr = VecGetArray(vec_from,&array_from);CHKERRQ(ierr); 554 } 555 if (vec_to) { 556 ierr = VecGetArray(vec_to,&array_to);CHKERRQ(ierr); 557 } 558 switch(pcbddc->coarse_problem_type){ 559 case SEQUENTIAL_BDDC: 560 if (smode == SCATTER_FORWARD) { 561 ierr = MPI_Gatherv(&array_from[0],pcbddc->local_primal_size,MPIU_SCALAR,&pcbddc->replicated_local_primal_values[0],pcbddc->local_primal_sizes,pcbddc->local_primal_displacements,MPIU_SCALAR,0,comm);CHKERRQ(ierr); 562 if (vec_to) { 563 if (imode == ADD_VALUES) { 564 for (i=0;i<pcbddc->replicated_primal_size;i++) { 565 array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i]; 566 } 567 } else { 568 for (i=0;i<pcbddc->replicated_primal_size;i++) { 569 array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i]; 570 } 571 } 572 } 573 } else { 574 if (vec_from) { 575 if (imode == ADD_VALUES) { 576 MPI_Comm vec_from_comm; 577 ierr = PetscObjectGetComm((PetscObject)(vec_from),&vec_from_comm);CHKERRQ(ierr); 578 SETERRQ2(vec_from_comm,PETSC_ERR_SUP,"Unsupported insert mode ADD_VALUES for SCATTER_REVERSE in %s for case %d\n",__FUNCT__,pcbddc->coarse_problem_type); 579 } 580 for (i=0;i<pcbddc->replicated_primal_size;i++) { 581 pcbddc->replicated_local_primal_values[i]=array_from[pcbddc->replicated_local_primal_indices[i]]; 582 } 583 } 584 ierr = MPI_Scatterv(&pcbddc->replicated_local_primal_values[0],pcbddc->local_primal_sizes,pcbddc->local_primal_displacements,MPIU_SCALAR,&array_to[0],pcbddc->local_primal_size,MPIU_SCALAR,0,comm);CHKERRQ(ierr); 585 } 586 break; 587 case REPLICATED_BDDC: 588 if (smode == SCATTER_FORWARD) { 589 ierr = MPI_Allgatherv(&array_from[0],pcbddc->local_primal_size,MPIU_SCALAR,&pcbddc->replicated_local_primal_values[0],pcbddc->local_primal_sizes,pcbddc->local_primal_displacements,MPIU_SCALAR,comm);CHKERRQ(ierr); 590 if (imode == ADD_VALUES) { 591 for (i=0;i<pcbddc->replicated_primal_size;i++) { 592 array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i]; 593 } 594 } else { 595 for (i=0;i<pcbddc->replicated_primal_size;i++) { 596 array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i]; 597 } 598 } 599 } else { /* no communications needed for SCATTER_REVERSE since needed data is already present */ 600 if (imode == ADD_VALUES) { 601 for (i=0;i<pcbddc->local_primal_size;i++) { 602 array_to[i]+=array_from[pcbddc->local_primal_indices[i]]; 603 } 604 } else { 605 for (i=0;i<pcbddc->local_primal_size;i++) { 606 array_to[i]=array_from[pcbddc->local_primal_indices[i]]; 607 } 608 } 609 } 610 break; 611 case MULTILEVEL_BDDC: 612 break; 613 case PARALLEL_BDDC: 614 break; 615 } 616 if (vec_from) { 617 ierr = VecRestoreArray(vec_from,&array_from);CHKERRQ(ierr); 618 } 619 if (vec_to) { 620 ierr = VecRestoreArray(vec_to,&array_to);CHKERRQ(ierr); 621 } 622 break; 623 } 624 PetscFunctionReturn(0); 625 } 626 627 #undef __FUNCT__ 628 #define __FUNCT__ "PCBDDCConstraintsSetUp" 629 PetscErrorCode PCBDDCConstraintsSetUp(PC pc) 630 { 631 PetscErrorCode ierr; 632 PC_IS* pcis = (PC_IS*)(pc->data); 633 PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 634 Mat_IS *matis = (Mat_IS*)pc->pmat->data; 635 PetscInt *nnz,*is_indices; 636 PetscScalar *temp_quadrature_constraint; 637 PetscInt *temp_indices,*temp_indices_to_constraint,*temp_indices_to_constraint_B,*local_to_B; 638 PetscInt local_primal_size,i,j,k,total_counts,max_size_of_constraint; 639 PetscInt n_vertices,size_of_constraint; 640 PetscReal real_value; 641 PetscBool nnsp_has_cnst=PETSC_FALSE,use_nnsp_true=pcbddc->use_nnsp_true; 642 PetscInt nnsp_size=0,nnsp_addone=0,temp_constraints,temp_start_ptr,n_ISForFaces,n_ISForEdges; 643 IS *used_IS,ISForVertices,*ISForFaces,*ISForEdges; 644 MatType impMatType=MATSEQAIJ; 645 PetscBLASInt Bs,Bt,lwork,lierr; 646 PetscReal tol=1.0e-8; 647 MatNullSpace nearnullsp; 648 const Vec *nearnullvecs; 649 Vec *localnearnullsp; 650 PetscScalar *work,*temp_basis,*array_vector,*correlation_mat; 651 PetscReal *rwork,*singular_vals; 652 PetscBLASInt Bone=1,*ipiv; 653 Vec temp_vec; 654 Mat temp_mat; 655 KSP temp_ksp; 656 PC temp_pc; 657 PetscInt s,start_constraint,dual_dofs; 658 PetscBool compute_submatrix,useksp=PETSC_FALSE; 659 PetscInt *aux_primal_permutation,*aux_primal_numbering; 660 PetscBool boolforchange,*change_basis; 661 /* some ugly conditional declarations */ 662 #if defined(PETSC_MISSING_LAPACK_GESVD) 663 PetscScalar one=1.0,zero=0.0; 664 PetscInt ii; 665 PetscScalar *singular_vectors; 666 PetscBLASInt *iwork,*ifail; 667 PetscReal dummy_real,abs_tol; 668 PetscBLASInt eigs_found; 669 #endif 670 PetscBLASInt dummy_int; 671 PetscScalar dummy_scalar; 672 PetscBool used_vertex,get_faces,get_edges,get_vertices; 673 674 PetscFunctionBegin; 675 /* Get index sets for faces, edges and vertices from graph */ 676 get_faces = PETSC_TRUE; 677 get_edges = PETSC_TRUE; 678 get_vertices = PETSC_TRUE; 679 if (pcbddc->vertices_flag) { 680 get_faces = PETSC_FALSE; 681 get_edges = PETSC_FALSE; 682 } 683 if (pcbddc->constraints_flag) { 684 get_vertices = PETSC_FALSE; 685 } 686 if (pcbddc->faces_flag) { 687 get_edges = PETSC_FALSE; 688 } 689 if (pcbddc->edges_flag) { 690 get_faces = PETSC_FALSE; 691 } 692 /* default */ 693 if (!get_faces && !get_edges && !get_vertices) { 694 get_vertices = PETSC_TRUE; 695 } 696 ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,get_faces,get_edges,get_vertices,&n_ISForFaces,&ISForFaces,&n_ISForEdges,&ISForEdges,&ISForVertices); 697 if (pcbddc->dbg_flag) { 698 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr); 699 i = 0; 700 if (ISForVertices) { 701 ierr = ISGetSize(ISForVertices,&i);CHKERRQ(ierr); 702 } 703 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate vertices\n",PetscGlobalRank,i);CHKERRQ(ierr); 704 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate edges\n",PetscGlobalRank,n_ISForEdges);CHKERRQ(ierr); 705 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate faces\n",PetscGlobalRank,n_ISForFaces);CHKERRQ(ierr); 706 ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 707 } 708 /* check if near null space is attached to global mat */ 709 ierr = MatGetNearNullSpace(pc->pmat,&nearnullsp);CHKERRQ(ierr); 710 if (nearnullsp) { 711 ierr = MatNullSpaceGetVecs(nearnullsp,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr); 712 } else { /* if near null space is not provided it uses constants */ 713 nnsp_has_cnst = PETSC_TRUE; 714 use_nnsp_true = PETSC_TRUE; 715 } 716 if (nnsp_has_cnst) { 717 nnsp_addone = 1; 718 } 719 /* 720 Evaluate maximum storage size needed by the procedure 721 - temp_indices will contain start index of each constraint stored as follows 722 - temp_indices_to_constraint [temp_indices[i],...,temp[indices[i+1]-1] will contain the indices (in local numbering) on which the constraint acts 723 - temp_indices_to_constraint_B[temp_indices[i],...,temp[indices[i+1]-1] will contain the indices (in boundary numbering) on which the constraint acts 724 - temp_quadrature_constraint [temp_indices[i],...,temp[indices[i+1]-1] will contain the scalars representing the constraint itself 725 */ 726 total_counts = n_ISForFaces+n_ISForEdges; 727 total_counts *= (nnsp_addone+nnsp_size); 728 n_vertices = 0; 729 if (ISForVertices) { 730 ierr = ISGetSize(ISForVertices,&n_vertices);CHKERRQ(ierr); 731 } 732 total_counts += n_vertices; 733 ierr = PetscMalloc((total_counts+1)*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr); 734 ierr = PetscMalloc((total_counts+1)*sizeof(PetscBool),&change_basis);CHKERRQ(ierr); 735 total_counts = 0; 736 max_size_of_constraint = 0; 737 for (i=0;i<n_ISForEdges+n_ISForFaces;i++) { 738 if (i<n_ISForEdges) { 739 used_IS = &ISForEdges[i]; 740 } else { 741 used_IS = &ISForFaces[i-n_ISForEdges]; 742 } 743 ierr = ISGetSize(*used_IS,&j);CHKERRQ(ierr); 744 total_counts += j; 745 max_size_of_constraint = PetscMax(j,max_size_of_constraint); 746 } 747 total_counts *= (nnsp_addone+nnsp_size); 748 total_counts += n_vertices; 749 ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&temp_quadrature_constraint);CHKERRQ(ierr); 750 ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint);CHKERRQ(ierr); 751 ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint_B);CHKERRQ(ierr); 752 ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&local_to_B);CHKERRQ(ierr); 753 ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 754 for (i=0;i<pcis->n;i++) { 755 local_to_B[i]=-1; 756 } 757 for (i=0;i<pcis->n_B;i++) { 758 local_to_B[is_indices[i]]=i; 759 } 760 ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 761 762 /* First we issue queries to allocate optimal workspace for LAPACKgesvd or LAPACKsyev/LAPACKheev */ 763 rwork = 0; 764 work = 0; 765 singular_vals = 0; 766 temp_basis = 0; 767 correlation_mat = 0; 768 if (!pcbddc->use_nnsp_true) { 769 PetscScalar temp_work; 770 #if defined(PETSC_MISSING_LAPACK_GESVD) 771 /* POD */ 772 PetscInt max_n; 773 max_n = nnsp_addone+nnsp_size; 774 /* using some techniques borrowed from Proper Orthogonal Decomposition */ 775 ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&correlation_mat);CHKERRQ(ierr); 776 ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&singular_vectors);CHKERRQ(ierr); 777 ierr = PetscMalloc(max_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr); 778 ierr = PetscMalloc(max_size_of_constraint*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr); 779 #if defined(PETSC_USE_COMPLEX) 780 ierr = PetscMalloc(3*max_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr); 781 #endif 782 ierr = PetscMalloc(5*max_n*sizeof(PetscBLASInt),&iwork);CHKERRQ(ierr); 783 ierr = PetscMalloc(max_n*sizeof(PetscBLASInt),&ifail);CHKERRQ(ierr); 784 /* now we evaluate the optimal workspace using query with lwork=-1 */ 785 ierr = PetscBLASIntCast(max_n,&Bt);CHKERRQ(ierr); 786 lwork=-1; 787 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 788 #if !defined(PETSC_USE_COMPLEX) 789 abs_tol=1.e-8; 790 /* LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,&lierr); */ 791 PetscStackCallBLAS("LAPACKsyevx",LAPACKsyevx_("V","A","U",&Bt,correlation_mat,&Bt,&dummy_real,&dummy_real,&dummy_int,&dummy_int,&abs_tol,&eigs_found,singular_vals,singular_vectors,&Bt,&temp_work,&lwork,iwork,ifail,&lierr)); 792 #else 793 /* LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,rwork,&lierr); */ 794 /* LAPACK call is missing here! TODO */ 795 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1"); 796 #endif 797 if ( lierr ) { 798 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYEVX Lapack routine %d",(int)lierr); 799 } 800 ierr = PetscFPTrapPop();CHKERRQ(ierr); 801 #else /* on missing GESVD */ 802 /* SVD */ 803 PetscInt max_n,min_n; 804 max_n = max_size_of_constraint; 805 min_n = nnsp_addone+nnsp_size; 806 if (max_size_of_constraint < ( nnsp_addone+nnsp_size ) ) { 807 min_n = max_size_of_constraint; 808 max_n = nnsp_addone+nnsp_size; 809 } 810 ierr = PetscMalloc(min_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr); 811 #if defined(PETSC_USE_COMPLEX) 812 ierr = PetscMalloc(5*min_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr); 813 #endif 814 /* now we evaluate the optimal workspace using query with lwork=-1 */ 815 lwork=-1; 816 ierr = PetscBLASIntCast(max_n,&Bs);CHKERRQ(ierr); 817 ierr = PetscBLASIntCast(min_n,&Bt);CHKERRQ(ierr); 818 dummy_int = Bs; 819 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 820 #if !defined(PETSC_USE_COMPLEX) 821 PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Bs,singular_vals,&dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,&temp_work,&lwork,&lierr)); 822 #else 823 PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Bs,singular_vals,&dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,&temp_work,&lwork,rwork,&lierr)); 824 #endif 825 if ( lierr ) { 826 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SVD Lapack routine %d",(int)lierr); 827 } 828 ierr = PetscFPTrapPop();CHKERRQ(ierr); 829 #endif 830 /* Allocate optimal workspace */ 831 ierr = PetscBLASIntCast((PetscInt)PetscRealPart(temp_work),&lwork);CHKERRQ(ierr); 832 total_counts = (PetscInt)lwork; 833 ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&work);CHKERRQ(ierr); 834 } 835 /* get local part of global near null space vectors */ 836 ierr = PetscMalloc(nnsp_size*sizeof(Vec),&localnearnullsp);CHKERRQ(ierr); 837 for (k=0;k<nnsp_size;k++) { 838 ierr = VecDuplicate(pcis->vec1_N,&localnearnullsp[k]);CHKERRQ(ierr); 839 ierr = VecScatterBegin(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 840 ierr = VecScatterEnd(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 841 } 842 /* Now we can loop on constraining sets */ 843 total_counts = 0; 844 temp_indices[0] = 0; 845 /* vertices */ 846 if (ISForVertices) { 847 ierr = ISGetIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); 848 if (nnsp_has_cnst) { /* consider all vertices */ 849 for (i=0;i<n_vertices;i++) { 850 temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i]; 851 temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]]; 852 temp_quadrature_constraint[temp_indices[total_counts]]=1.0; 853 temp_indices[total_counts+1]=temp_indices[total_counts]+1; 854 change_basis[total_counts]=PETSC_FALSE; 855 total_counts++; 856 } 857 } else { /* consider vertices for which exist at least a localnearnullsp which is not null there */ 858 for (i=0;i<n_vertices;i++) { 859 used_vertex=PETSC_FALSE; 860 k=0; 861 while (!used_vertex && k<nnsp_size) { 862 ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 863 if (PetscAbsScalar(array_vector[is_indices[i]])>0.0) { 864 temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i]; 865 temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]]; 866 temp_quadrature_constraint[temp_indices[total_counts]]=1.0; 867 temp_indices[total_counts+1]=temp_indices[total_counts]+1; 868 change_basis[total_counts]=PETSC_FALSE; 869 total_counts++; 870 used_vertex=PETSC_TRUE; 871 } 872 ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 873 k++; 874 } 875 } 876 } 877 ierr = ISRestoreIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); 878 n_vertices = total_counts; 879 } 880 /* edges and faces */ 881 for (i=0;i<n_ISForEdges+n_ISForFaces;i++) { 882 if (i<n_ISForEdges) { 883 used_IS = &ISForEdges[i]; 884 boolforchange = pcbddc->use_change_of_basis; 885 } else { 886 used_IS = &ISForFaces[i-n_ISForEdges]; 887 boolforchange = pcbddc->use_change_on_faces; 888 } 889 temp_constraints = 0; /* zero the number of constraints I have on this conn comp */ 890 temp_start_ptr = total_counts; /* need to know the starting index of constraints stored */ 891 ierr = ISGetSize(*used_IS,&size_of_constraint);CHKERRQ(ierr); 892 ierr = ISGetIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 893 /* HACK: change of basis should not performed on local periodic nodes */ 894 if (pcbddc->mat_graph->mirrors && pcbddc->mat_graph->mirrors[is_indices[0]]) { 895 boolforchange = PETSC_FALSE; 896 } 897 if (nnsp_has_cnst) { 898 PetscScalar quad_value; 899 temp_constraints++; 900 quad_value = (PetscScalar)(1.0/PetscSqrtReal((PetscReal)size_of_constraint)); 901 for (j=0;j<size_of_constraint;j++) { 902 temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j]; 903 temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]]; 904 temp_quadrature_constraint[temp_indices[total_counts]+j]=quad_value; 905 } 906 temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint; /* store new starting point */ 907 change_basis[total_counts]=boolforchange; 908 total_counts++; 909 } 910 for (k=0;k<nnsp_size;k++) { 911 ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 912 for (j=0;j<size_of_constraint;j++) { 913 temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j]; 914 temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]]; 915 temp_quadrature_constraint[temp_indices[total_counts]+j]=array_vector[is_indices[j]]; 916 } 917 ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 918 real_value = 1.0; 919 if (use_nnsp_true) { /* check if array is null on the connected component in case use_nnsp_true has been requested */ 920 ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr); 921 PetscStackCallBLAS("BLASasum",real_value = BLASasum_(&Bs,&temp_quadrature_constraint[temp_indices[total_counts]],&Bone)); 922 } 923 if (real_value > 0.0) { /* keep indices and values */ 924 temp_constraints++; 925 temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint; /* store new starting point */ 926 change_basis[total_counts]=boolforchange; 927 total_counts++; 928 } 929 } 930 ierr = ISRestoreIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 931 /* perform SVD on the constraint if use_nnsp_true has not be requested by the user */ 932 if (!use_nnsp_true) { 933 ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr); 934 ierr = PetscBLASIntCast(temp_constraints,&Bt);CHKERRQ(ierr); 935 936 #if defined(PETSC_MISSING_LAPACK_GESVD) 937 ierr = PetscMemzero(correlation_mat,Bt*Bt*sizeof(PetscScalar));CHKERRQ(ierr); 938 /* Store upper triangular part of correlation matrix */ 939 for (j=0;j<temp_constraints;j++) { 940 for (k=0;k<j+1;k++) { 941 PetscStackCallBLAS("BLASdot",correlation_mat[j*temp_constraints+k]=BLASdot_(&Bs,&temp_quadrature_constraint[temp_indices[temp_start_ptr+j]],&Bone,&temp_quadrature_constraint[temp_indices[temp_start_ptr+k]],&Bone)); 942 943 } 944 } 945 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 946 #if !defined(PETSC_USE_COMPLEX) 947 /* LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,work,&lwork,&lierr); */ 948 PetscStackCallBLAS("LAPACKsyevx",LAPACKsyevx_("V","A","U",&Bt,correlation_mat,&Bt,&dummy_real,&dummy_real,&dummy_int,&dummy_int,&abs_tol,&eigs_found,singular_vals,singular_vectors,&Bt,work,&lwork,iwork,ifail,&lierr)); 949 #else 950 /* LAPACK call is missing here! TODO */ 951 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1"); 952 #endif 953 if (lierr) { 954 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYEVX Lapack routine %d",(int)lierr); 955 } 956 ierr = PetscFPTrapPop();CHKERRQ(ierr); 957 /* retain eigenvalues greater than tol: note that lapack SYEV gives eigs in ascending order */ 958 j=0; 959 while (j < Bt && singular_vals[j] < tol) j++; 960 total_counts=total_counts-j; 961 if (j<temp_constraints) { 962 for (k=j;k<Bt;k++) { 963 singular_vals[k]=1.0/PetscSqrtReal(singular_vals[k]); 964 } 965 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 966 PetscStackCallBLAS("BLASgemm_",BLASgemm_("N","N",&Bs,&Bt,&Bt,&one,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,correlation_mat,&Bt,&zero,temp_basis,&Bs)); 967 ierr = PetscFPTrapPop();CHKERRQ(ierr); 968 /* copy POD basis into used quadrature memory */ 969 for (k=0;k<Bt-j;k++) { 970 for (ii=0;ii<size_of_constraint;ii++) { 971 temp_quadrature_constraint[temp_indices[temp_start_ptr+k]+ii]=singular_vals[Bt-1-k]*temp_basis[(Bt-1-k)*size_of_constraint+ii]; 972 } 973 } 974 } 975 976 #else /* on missing GESVD */ 977 PetscInt min_n = temp_constraints; 978 if (min_n > size_of_constraint) min_n = size_of_constraint; 979 dummy_int = Bs; 980 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 981 #if !defined(PETSC_USE_COMPLEX) 982 PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,singular_vals,&dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,work,&lwork,&lierr)); 983 #else 984 PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,singular_vals,&dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,work,&lwork,rwork,&lierr)); 985 #endif 986 if (lierr) { 987 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SVD Lapack routine %d",(int)lierr); 988 } 989 ierr = PetscFPTrapPop();CHKERRQ(ierr); 990 /* retain eigenvalues greater than tol: note that lapack SVD gives eigs in descending order */ 991 j = 0; 992 while (j < min_n && singular_vals[min_n-j-1] < tol) j++; 993 total_counts = total_counts-(PetscInt)Bt+(min_n-j); 994 #endif 995 } 996 } 997 /* free index sets of faces, edges and vertices */ 998 for (i=0;i<n_ISForFaces;i++) { 999 ierr = ISDestroy(&ISForFaces[i]);CHKERRQ(ierr); 1000 } 1001 ierr = PetscFree(ISForFaces);CHKERRQ(ierr); 1002 for (i=0;i<n_ISForEdges;i++) { 1003 ierr = ISDestroy(&ISForEdges[i]);CHKERRQ(ierr); 1004 } 1005 ierr = PetscFree(ISForEdges);CHKERRQ(ierr); 1006 ierr = ISDestroy(&ISForVertices);CHKERRQ(ierr); 1007 1008 /* set quantities in pcbddc data structure */ 1009 /* n_vertices defines the number of point primal dofs */ 1010 /* n_constraints defines the number of averages (they can be point primal dofs if change of basis is requested) */ 1011 local_primal_size = total_counts; 1012 pcbddc->n_vertices = n_vertices; 1013 pcbddc->n_constraints = total_counts-n_vertices; 1014 pcbddc->local_primal_size = local_primal_size; 1015 1016 /* Create constraint matrix */ 1017 /* The constraint matrix is used to compute the l2g map of primal dofs */ 1018 /* so we need to set it up properly either with or without change of basis */ 1019 ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ConstraintMatrix);CHKERRQ(ierr); 1020 ierr = MatSetType(pcbddc->ConstraintMatrix,impMatType);CHKERRQ(ierr); 1021 ierr = MatSetSizes(pcbddc->ConstraintMatrix,local_primal_size,pcis->n,local_primal_size,pcis->n);CHKERRQ(ierr); 1022 /* compute a local numbering of constraints : vertices first then constraints */ 1023 ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr); 1024 ierr = VecGetArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr); 1025 ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_numbering);CHKERRQ(ierr); 1026 ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_permutation);CHKERRQ(ierr); 1027 total_counts=0; 1028 /* find vertices: subdomain corners plus dofs with basis changed */ 1029 for (i=0;i<local_primal_size;i++) { 1030 size_of_constraint=temp_indices[i+1]-temp_indices[i]; 1031 if (change_basis[i] || size_of_constraint == 1) { 1032 k=0; 1033 while(k < size_of_constraint && array_vector[temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1]] != 0.0) { 1034 k=k+1; 1035 } 1036 j=temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1]; 1037 array_vector[j] = 1.0; 1038 aux_primal_numbering[total_counts]=j; 1039 aux_primal_permutation[total_counts]=total_counts; 1040 total_counts++; 1041 } 1042 } 1043 ierr = VecRestoreArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr); 1044 /* permute indices in order to have a sorted set of vertices */ 1045 ierr = PetscSortIntWithPermutation(total_counts,aux_primal_numbering,aux_primal_permutation); 1046 /* nonzero structure */ 1047 ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 1048 for (i=0;i<total_counts;i++) { 1049 nnz[i]=1; 1050 } 1051 j=total_counts; 1052 for (i=n_vertices;i<local_primal_size;i++) { 1053 if (!change_basis[i]) { 1054 nnz[j]=temp_indices[i+1]-temp_indices[i]; 1055 j++; 1056 } 1057 } 1058 ierr = MatSeqAIJSetPreallocation(pcbddc->ConstraintMatrix,0,nnz);CHKERRQ(ierr); 1059 ierr = PetscFree(nnz);CHKERRQ(ierr); 1060 /* set values in constraint matrix */ 1061 for (i=0;i<total_counts;i++) { 1062 j = aux_primal_permutation[i]; 1063 k = aux_primal_numbering[j]; 1064 ierr = MatSetValue(pcbddc->ConstraintMatrix,i,k,1.0,INSERT_VALUES);CHKERRQ(ierr); 1065 } 1066 for (i=n_vertices;i<local_primal_size;i++) { 1067 if (!change_basis[i]) { 1068 size_of_constraint=temp_indices[i+1]-temp_indices[i]; 1069 ierr = MatSetValues(pcbddc->ConstraintMatrix,1,&total_counts,size_of_constraint,&temp_indices_to_constraint[temp_indices[i]],&temp_quadrature_constraint[temp_indices[i]],INSERT_VALUES);CHKERRQ(ierr); 1070 total_counts++; 1071 } 1072 } 1073 ierr = PetscFree(aux_primal_numbering);CHKERRQ(ierr); 1074 ierr = PetscFree(aux_primal_permutation);CHKERRQ(ierr); 1075 /* assembling */ 1076 ierr = MatAssemblyBegin(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1077 ierr = MatAssemblyEnd(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1078 1079 /* Create matrix for change of basis. We don't need it in case pcbddc->use_change_of_basis is FALSE */ 1080 if (pcbddc->use_change_of_basis) { 1081 ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr); 1082 ierr = MatSetType(pcbddc->ChangeOfBasisMatrix,impMatType);CHKERRQ(ierr); 1083 ierr = MatSetSizes(pcbddc->ChangeOfBasisMatrix,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);CHKERRQ(ierr); 1084 /* work arrays */ 1085 /* we need to reuse these arrays, so we free them */ 1086 ierr = PetscFree(temp_basis);CHKERRQ(ierr); 1087 ierr = PetscFree(work);CHKERRQ(ierr); 1088 ierr = PetscMalloc(pcis->n_B*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 1089 ierr = PetscMalloc((nnsp_addone+nnsp_size)*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr); 1090 ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscScalar),&work);CHKERRQ(ierr); 1091 ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscBLASInt),&ipiv);CHKERRQ(ierr); 1092 for (i=0;i<pcis->n_B;i++) { 1093 nnz[i]=1; 1094 } 1095 /* Overestimated nonzeros per row */ 1096 k=1; 1097 for (i=pcbddc->n_vertices;i<local_primal_size;i++) { 1098 if (change_basis[i]) { 1099 size_of_constraint = temp_indices[i+1]-temp_indices[i]; 1100 if (k < size_of_constraint) { 1101 k = size_of_constraint; 1102 } 1103 for (j=0;j<size_of_constraint;j++) { 1104 nnz[temp_indices_to_constraint_B[temp_indices[i]+j]] = size_of_constraint; 1105 } 1106 } 1107 } 1108 ierr = MatSeqAIJSetPreallocation(pcbddc->ChangeOfBasisMatrix,0,nnz);CHKERRQ(ierr); 1109 ierr = PetscFree(nnz);CHKERRQ(ierr); 1110 /* Temporary array to store indices */ 1111 ierr = PetscMalloc(k*sizeof(PetscInt),&is_indices);CHKERRQ(ierr); 1112 /* Set initial identity in the matrix */ 1113 for (i=0;i<pcis->n_B;i++) { 1114 ierr = MatSetValue(pcbddc->ChangeOfBasisMatrix,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr); 1115 } 1116 /* Now we loop on the constraints which need a change of basis */ 1117 /* Change of basis matrix is evaluated as the FIRST APPROACH in */ 1118 /* Klawonn and Widlund, Dual-primal FETI-DP methods for linear elasticity, (6.2.1) */ 1119 temp_constraints = 0; 1120 if (pcbddc->n_vertices < local_primal_size) { 1121 temp_start_ptr = temp_indices_to_constraint_B[temp_indices[pcbddc->n_vertices]]; 1122 } 1123 for (i=pcbddc->n_vertices;i<local_primal_size;i++) { 1124 if (change_basis[i]) { 1125 compute_submatrix = PETSC_FALSE; 1126 useksp = PETSC_FALSE; 1127 if (temp_start_ptr == temp_indices_to_constraint_B[temp_indices[i]]) { 1128 temp_constraints++; 1129 if (i == local_primal_size -1 || temp_start_ptr != temp_indices_to_constraint_B[temp_indices[i+1]]) { 1130 compute_submatrix = PETSC_TRUE; 1131 } 1132 } 1133 if (compute_submatrix) { 1134 if (temp_constraints > 1 || pcbddc->use_nnsp_true) { 1135 useksp = PETSC_TRUE; 1136 } 1137 size_of_constraint = temp_indices[i+1]-temp_indices[i]; 1138 if (useksp) { /* experimental TODO: reuse KSP and MAT instead of creating them each time */ 1139 ierr = MatCreate(PETSC_COMM_SELF,&temp_mat);CHKERRQ(ierr); 1140 ierr = MatSetType(temp_mat,impMatType);CHKERRQ(ierr); 1141 ierr = MatSetSizes(temp_mat,size_of_constraint,size_of_constraint,size_of_constraint,size_of_constraint);CHKERRQ(ierr); 1142 ierr = MatSeqAIJSetPreallocation(temp_mat,size_of_constraint,NULL);CHKERRQ(ierr); 1143 } 1144 /* First _size_of_constraint-temp_constraints_ columns */ 1145 dual_dofs = size_of_constraint-temp_constraints; 1146 start_constraint = i+1-temp_constraints; 1147 for (s=0;s<dual_dofs;s++) { 1148 is_indices[0] = s; 1149 for (j=0;j<temp_constraints;j++) { 1150 for (k=0;k<temp_constraints;k++) { 1151 temp_basis[j*temp_constraints+k]=temp_quadrature_constraint[temp_indices[start_constraint+k]+s+j+1]; 1152 } 1153 work[j]=-temp_quadrature_constraint[temp_indices[start_constraint+j]+s]; 1154 is_indices[j+1]=s+j+1; 1155 } 1156 Bt = temp_constraints; 1157 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1158 PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&Bt,&Bone,temp_basis,&Bt,ipiv,work,&Bt,&lierr)); 1159 if ( lierr ) { 1160 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GESV Lapack routine %d",(int)lierr); 1161 } 1162 ierr = PetscFPTrapPop();CHKERRQ(ierr); 1163 j = temp_indices_to_constraint_B[temp_indices[start_constraint]+s]; 1164 ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,temp_constraints,&temp_indices_to_constraint_B[temp_indices[start_constraint]+s+1],1,&j,work,INSERT_VALUES);CHKERRQ(ierr); 1165 if (useksp) { 1166 /* temp mat with transposed rows and columns */ 1167 ierr = MatSetValues(temp_mat,1,&s,temp_constraints,&is_indices[1],work,INSERT_VALUES);CHKERRQ(ierr); 1168 ierr = MatSetValue(temp_mat,is_indices[0],is_indices[0],1.0,INSERT_VALUES);CHKERRQ(ierr); 1169 } 1170 } 1171 if (useksp) { 1172 /* last rows of temp_mat */ 1173 for (j=0;j<size_of_constraint;j++) { 1174 is_indices[j] = j; 1175 } 1176 for (s=0;s<temp_constraints;s++) { 1177 k = s + dual_dofs; 1178 ierr = MatSetValues(temp_mat,1,&k,size_of_constraint,is_indices,&temp_quadrature_constraint[temp_indices[start_constraint+s]],INSERT_VALUES);CHKERRQ(ierr); 1179 } 1180 ierr = MatAssemblyBegin(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1181 ierr = MatAssemblyEnd(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1182 ierr = MatGetVecs(temp_mat,&temp_vec,NULL);CHKERRQ(ierr); 1183 ierr = KSPCreate(PETSC_COMM_SELF,&temp_ksp);CHKERRQ(ierr); 1184 ierr = KSPSetOperators(temp_ksp,temp_mat,temp_mat,SAME_PRECONDITIONER);CHKERRQ(ierr); 1185 ierr = KSPSetType(temp_ksp,KSPPREONLY);CHKERRQ(ierr); 1186 ierr = KSPGetPC(temp_ksp,&temp_pc);CHKERRQ(ierr); 1187 ierr = PCSetType(temp_pc,PCLU);CHKERRQ(ierr); 1188 ierr = KSPSetUp(temp_ksp);CHKERRQ(ierr); 1189 for (s=0;s<temp_constraints;s++) { 1190 ierr = VecSet(temp_vec,0.0);CHKERRQ(ierr); 1191 ierr = VecSetValue(temp_vec,s+dual_dofs,1.0,INSERT_VALUES);CHKERRQ(ierr); 1192 ierr = VecAssemblyBegin(temp_vec);CHKERRQ(ierr); 1193 ierr = VecAssemblyEnd(temp_vec);CHKERRQ(ierr); 1194 ierr = KSPSolve(temp_ksp,temp_vec,temp_vec);CHKERRQ(ierr); 1195 ierr = VecGetArray(temp_vec,&array_vector);CHKERRQ(ierr); 1196 j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1]; 1197 /* last columns of change of basis matrix associated to new primal dofs */ 1198 ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,&temp_indices_to_constraint_B[temp_indices[start_constraint+s]],1,&j,array_vector,INSERT_VALUES);CHKERRQ(ierr); 1199 ierr = VecRestoreArray(temp_vec,&array_vector);CHKERRQ(ierr); 1200 } 1201 ierr = MatDestroy(&temp_mat);CHKERRQ(ierr); 1202 ierr = KSPDestroy(&temp_ksp);CHKERRQ(ierr); 1203 ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 1204 } else { 1205 /* last columns of change of basis matrix associated to new primal dofs */ 1206 for (s=0;s<temp_constraints;s++) { 1207 j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1]; 1208 ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,&temp_indices_to_constraint_B[temp_indices[start_constraint+s]],1,&j,&temp_quadrature_constraint[temp_indices[start_constraint+s]],INSERT_VALUES);CHKERRQ(ierr); 1209 } 1210 } 1211 /* prepare for the next cycle */ 1212 temp_constraints = 0; 1213 if (i != local_primal_size -1 ) { 1214 temp_start_ptr = temp_indices_to_constraint_B[temp_indices[i+1]]; 1215 } 1216 } 1217 } 1218 } 1219 /* assembling */ 1220 ierr = MatAssemblyBegin(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1221 ierr = MatAssemblyEnd(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1222 ierr = PetscFree(ipiv);CHKERRQ(ierr); 1223 ierr = PetscFree(is_indices);CHKERRQ(ierr); 1224 } 1225 /* free workspace no longer needed */ 1226 ierr = PetscFree(rwork);CHKERRQ(ierr); 1227 ierr = PetscFree(work);CHKERRQ(ierr); 1228 ierr = PetscFree(temp_basis);CHKERRQ(ierr); 1229 ierr = PetscFree(singular_vals);CHKERRQ(ierr); 1230 ierr = PetscFree(correlation_mat);CHKERRQ(ierr); 1231 ierr = PetscFree(temp_indices);CHKERRQ(ierr); 1232 ierr = PetscFree(change_basis);CHKERRQ(ierr); 1233 ierr = PetscFree(temp_indices_to_constraint);CHKERRQ(ierr); 1234 ierr = PetscFree(temp_indices_to_constraint_B);CHKERRQ(ierr); 1235 ierr = PetscFree(local_to_B);CHKERRQ(ierr); 1236 ierr = PetscFree(temp_quadrature_constraint);CHKERRQ(ierr); 1237 #if defined(PETSC_MISSING_LAPACK_GESVD) 1238 ierr = PetscFree(iwork);CHKERRQ(ierr); 1239 ierr = PetscFree(ifail);CHKERRQ(ierr); 1240 ierr = PetscFree(singular_vectors);CHKERRQ(ierr); 1241 #endif 1242 for (k=0;k<nnsp_size;k++) { 1243 ierr = VecDestroy(&localnearnullsp[k]);CHKERRQ(ierr); 1244 } 1245 ierr = PetscFree(localnearnullsp);CHKERRQ(ierr); 1246 PetscFunctionReturn(0); 1247 } 1248 1249 #undef __FUNCT__ 1250 #define __FUNCT__ "PCBDDCAnalyzeInterface" 1251 PetscErrorCode PCBDDCAnalyzeInterface(PC pc) 1252 { 1253 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1254 PC_IS *pcis = (PC_IS*)pc->data; 1255 Mat_IS *matis = (Mat_IS*)pc->pmat->data; 1256 PetscInt bs,ierr,i,vertex_size; 1257 PetscViewer viewer=pcbddc->dbg_viewer; 1258 1259 PetscFunctionBegin; 1260 /* Init local Graph struct */ 1261 ierr = PCBDDCGraphInit(pcbddc->mat_graph,matis->mapping);CHKERRQ(ierr); 1262 1263 /* Check validity of the csr graph passed in by the user */ 1264 if (pcbddc->mat_graph->nvtxs_csr != pcbddc->mat_graph->nvtxs) { 1265 ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr); 1266 } 1267 /* Set default CSR adjacency of local dofs if not provided by the user with PCBDDCSetLocalAdjacencyGraph */ 1268 if (!pcbddc->mat_graph->xadj || !pcbddc->mat_graph->adjncy) { 1269 Mat mat_adj; 1270 const PetscInt *xadj,*adjncy; 1271 PetscBool flg_row=PETSC_TRUE; 1272 1273 ierr = MatConvert(matis->A,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr); 1274 ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr); 1275 if (!flg_row) { 1276 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__); 1277 } 1278 ierr = PCBDDCSetLocalAdjacencyGraph(pc,i,xadj,adjncy,PETSC_COPY_VALUES);CHKERRQ(ierr); 1279 ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr); 1280 if (!flg_row) { 1281 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 1282 } 1283 ierr = MatDestroy(&mat_adj);CHKERRQ(ierr); 1284 } 1285 1286 /* Set default dofs' splitting if no information has been provided by the user with PCBDDCSetDofsSplitting */ 1287 vertex_size = 1; 1288 if (!pcbddc->n_ISForDofs) { 1289 IS *custom_ISForDofs; 1290 1291 ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 1292 ierr = PetscMalloc(bs*sizeof(IS),&custom_ISForDofs);CHKERRQ(ierr); 1293 for (i=0;i<bs;i++) { 1294 ierr = ISCreateStride(PETSC_COMM_SELF,pcis->n/bs,i,bs,&custom_ISForDofs[i]);CHKERRQ(ierr); 1295 } 1296 ierr = PCBDDCSetDofsSplitting(pc,bs,custom_ISForDofs);CHKERRQ(ierr); 1297 /* remove my references to IS objects */ 1298 for (i=0;i<bs;i++) { 1299 ierr = ISDestroy(&custom_ISForDofs[i]);CHKERRQ(ierr); 1300 } 1301 ierr = PetscFree(custom_ISForDofs);CHKERRQ(ierr); 1302 } else { /* mat block size as vertex size (used for elasticity) */ 1303 ierr = MatGetBlockSize(matis->A,&vertex_size);CHKERRQ(ierr); 1304 } 1305 1306 /* Setup of Graph */ 1307 ierr = PCBDDCGraphSetUp(pcbddc->mat_graph,vertex_size,pcbddc->NeumannBoundaries,pcbddc->DirichletBoundaries,pcbddc->n_ISForDofs,pcbddc->ISForDofs,pcbddc->user_primal_vertices); 1308 1309 /* Graph's connected components analysis */ 1310 ierr = PCBDDCGraphComputeConnectedComponents(pcbddc->mat_graph);CHKERRQ(ierr); 1311 1312 /* print some info to stdout */ 1313 if (pcbddc->dbg_flag) { 1314 ierr = PCBDDCGraphASCIIView(pcbddc->mat_graph,pcbddc->dbg_flag,viewer); 1315 } 1316 PetscFunctionReturn(0); 1317 } 1318 1319 #undef __FUNCT__ 1320 #define __FUNCT__ "PCBDDCGetPrimalVerticesLocalIdx" 1321 PetscErrorCode PCBDDCGetPrimalVerticesLocalIdx(PC pc, PetscInt *n_vertices, PetscInt *vertices_idx[]) 1322 { 1323 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 1324 PetscInt *vertices,*row_cmat_indices,n,i,size_of_constraint,local_primal_size; 1325 PetscErrorCode ierr; 1326 1327 PetscFunctionBegin; 1328 n = 0; 1329 vertices = 0; 1330 if (pcbddc->ConstraintMatrix) { 1331 ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&i);CHKERRQ(ierr); 1332 for (i=0;i<local_primal_size;i++) { 1333 ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 1334 if (size_of_constraint == 1) n++; 1335 ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 1336 } 1337 if (vertices_idx) { 1338 ierr = PetscMalloc(n*sizeof(PetscInt),&vertices);CHKERRQ(ierr); 1339 n = 0; 1340 for (i=0;i<local_primal_size;i++) { 1341 ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1342 if (size_of_constraint == 1) { 1343 vertices[n++]=row_cmat_indices[0]; 1344 } 1345 ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1346 } 1347 } 1348 } 1349 *n_vertices = n; 1350 if (vertices_idx) *vertices_idx = vertices; 1351 PetscFunctionReturn(0); 1352 } 1353 1354 #undef __FUNCT__ 1355 #define __FUNCT__ "PCBDDCGetPrimalConstraintsLocalIdx" 1356 PetscErrorCode PCBDDCGetPrimalConstraintsLocalIdx(PC pc, PetscInt *n_constraints, PetscInt *constraints_idx[]) 1357 { 1358 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 1359 PetscInt *constraints_index,*row_cmat_indices,*row_cmat_global_indices; 1360 PetscInt n,i,j,size_of_constraint,local_primal_size,local_size,max_size_of_constraint,min_index,min_loc; 1361 PetscBool *touched; 1362 PetscErrorCode ierr; 1363 1364 PetscFunctionBegin; 1365 n = 0; 1366 constraints_index = 0; 1367 if (pcbddc->ConstraintMatrix) { 1368 ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&local_size);CHKERRQ(ierr); 1369 max_size_of_constraint = 0; 1370 for (i=0;i<local_primal_size;i++) { 1371 ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 1372 if (size_of_constraint > 1) { 1373 n++; 1374 } 1375 max_size_of_constraint = PetscMax(size_of_constraint,max_size_of_constraint); 1376 ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 1377 } 1378 if (constraints_idx) { 1379 ierr = PetscMalloc(n*sizeof(PetscInt),&constraints_index);CHKERRQ(ierr); 1380 ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&row_cmat_global_indices);CHKERRQ(ierr); 1381 ierr = PetscMalloc(local_size*sizeof(PetscBool),&touched);CHKERRQ(ierr); 1382 ierr = PetscMemzero(touched,local_size*sizeof(PetscBool));CHKERRQ(ierr); 1383 n = 0; 1384 for (i=0;i<local_primal_size;i++) { 1385 ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1386 if (size_of_constraint > 1) { 1387 ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,row_cmat_indices,row_cmat_global_indices);CHKERRQ(ierr); 1388 /* find first untouched local node */ 1389 j = 0; 1390 while(touched[row_cmat_indices[j]]) j++; 1391 min_index = row_cmat_global_indices[j]; 1392 min_loc = j; 1393 /* search the minimum among nodes not yet touched on the connected component 1394 since there can be more than one constraint on a single cc */ 1395 for (j=1;j<size_of_constraint;j++) { 1396 if (min_index > row_cmat_global_indices[j] && !touched[row_cmat_indices[j]]) { 1397 min_index = row_cmat_global_indices[j]; 1398 min_loc = j; 1399 } 1400 } 1401 touched[row_cmat_indices[min_loc]] = PETSC_TRUE; 1402 constraints_index[n++] = row_cmat_indices[min_loc]; 1403 } 1404 ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1405 } 1406 ierr = PetscFree(touched);CHKERRQ(ierr); 1407 ierr = PetscFree(row_cmat_global_indices);CHKERRQ(ierr); 1408 } 1409 } 1410 *n_constraints = n; 1411 if (constraints_idx) *constraints_idx = constraints_index; 1412 PetscFunctionReturn(0); 1413 } 1414 1415 /* the next two functions has been adapted from pcis.c */ 1416 #undef __FUNCT__ 1417 #define __FUNCT__ "PCBDDCApplySchur" 1418 PetscErrorCode PCBDDCApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 1419 { 1420 PetscErrorCode ierr; 1421 PC_IS *pcis = (PC_IS*)(pc->data); 1422 1423 PetscFunctionBegin; 1424 if (!vec2_B) { vec2_B = v; } 1425 ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 1426 ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr); 1427 ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 1428 ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr); 1429 ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 1430 PetscFunctionReturn(0); 1431 } 1432 1433 #undef __FUNCT__ 1434 #define __FUNCT__ "PCBDDCApplySchurTranspose" 1435 PetscErrorCode PCBDDCApplySchurTranspose(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 1436 { 1437 PetscErrorCode ierr; 1438 PC_IS *pcis = (PC_IS*)(pc->data); 1439 1440 PetscFunctionBegin; 1441 if (!vec2_B) { vec2_B = v; } 1442 ierr = MatMultTranspose(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 1443 ierr = MatMultTranspose(pcis->A_BI,v,vec1_D);CHKERRQ(ierr); 1444 ierr = KSPSolveTranspose(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 1445 ierr = MatMultTranspose(pcis->A_IB,vec2_D,vec2_B);CHKERRQ(ierr); 1446 ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 1447 PetscFunctionReturn(0); 1448 } 1449 1450 #undef __FUNCT__ 1451 #define __FUNCT__ "PCBDDCSubsetNumbering" 1452 PetscErrorCode PCBDDCSubsetNumbering(MPI_Comm comm,ISLocalToGlobalMapping l2gmap, PetscInt n_local_dofs, PetscInt local_dofs[], PetscInt local_dofs_mult[], PetscInt* n_global_subset, PetscInt* global_numbering_subset[]) 1453 { 1454 Vec local_vec,global_vec; 1455 IS seqis,paris; 1456 VecScatter scatter_ctx; 1457 PetscScalar *array; 1458 PetscInt *temp_global_dofs; 1459 PetscScalar globalsum; 1460 PetscInt i,j,s; 1461 PetscInt nlocals,first_index,old_index,max_local; 1462 PetscMPIInt rank_prec_comm,size_prec_comm,max_global; 1463 PetscMPIInt *dof_sizes,*dof_displs; 1464 PetscBool first_found; 1465 PetscErrorCode ierr; 1466 1467 PetscFunctionBegin; 1468 /* mpi buffers */ 1469 MPI_Comm_size(comm,&size_prec_comm); 1470 MPI_Comm_rank(comm,&rank_prec_comm); 1471 j = ( !rank_prec_comm ? size_prec_comm : 0); 1472 ierr = PetscMalloc(j*sizeof(*dof_sizes),&dof_sizes);CHKERRQ(ierr); 1473 ierr = PetscMalloc(j*sizeof(*dof_displs),&dof_displs);CHKERRQ(ierr); 1474 /* get maximum size of subset */ 1475 ierr = PetscMalloc(n_local_dofs*sizeof(PetscInt),&temp_global_dofs);CHKERRQ(ierr); 1476 ierr = ISLocalToGlobalMappingApply(l2gmap,n_local_dofs,local_dofs,temp_global_dofs);CHKERRQ(ierr); 1477 max_local = 0; 1478 if (n_local_dofs) { 1479 max_local = temp_global_dofs[0]; 1480 for (i=1;i<n_local_dofs;i++) { 1481 if (max_local < temp_global_dofs[i] ) { 1482 max_local = temp_global_dofs[i]; 1483 } 1484 } 1485 } 1486 ierr = MPI_Allreduce(&max_local,&max_global,1,MPIU_INT,MPI_MAX,comm); 1487 max_global++; 1488 max_local = 0; 1489 if (n_local_dofs) { 1490 max_local = local_dofs[0]; 1491 for (i=1;i<n_local_dofs;i++) { 1492 if (max_local < local_dofs[i] ) { 1493 max_local = local_dofs[i]; 1494 } 1495 } 1496 } 1497 max_local++; 1498 /* allocate workspace */ 1499 ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr); 1500 ierr = VecSetSizes(local_vec,PETSC_DECIDE,max_local);CHKERRQ(ierr); 1501 ierr = VecSetType(local_vec,VECSEQ);CHKERRQ(ierr); 1502 ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr); 1503 ierr = VecSetSizes(global_vec,PETSC_DECIDE,max_global);CHKERRQ(ierr); 1504 ierr = VecSetType(global_vec,VECMPI);CHKERRQ(ierr); 1505 /* create scatter */ 1506 ierr = ISCreateGeneral(PETSC_COMM_SELF,n_local_dofs,local_dofs,PETSC_COPY_VALUES,&seqis);CHKERRQ(ierr); 1507 ierr = ISCreateGeneral(comm,n_local_dofs,temp_global_dofs,PETSC_COPY_VALUES,&paris);CHKERRQ(ierr); 1508 ierr = VecScatterCreate(local_vec,seqis,global_vec,paris,&scatter_ctx);CHKERRQ(ierr); 1509 ierr = ISDestroy(&seqis);CHKERRQ(ierr); 1510 ierr = ISDestroy(&paris);CHKERRQ(ierr); 1511 /* init array */ 1512 ierr = VecSet(global_vec,0.0);CHKERRQ(ierr); 1513 ierr = VecSet(local_vec,0.0);CHKERRQ(ierr); 1514 ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 1515 if (local_dofs_mult) { 1516 for (i=0;i<n_local_dofs;i++) { 1517 array[local_dofs[i]]=(PetscScalar)local_dofs_mult[i]; 1518 } 1519 } else { 1520 for (i=0;i<n_local_dofs;i++) { 1521 array[local_dofs[i]]=1.0; 1522 } 1523 } 1524 ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 1525 /* scatter into global vec and get total number of global dofs */ 1526 ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1527 ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1528 ierr = VecSum(global_vec,&globalsum);CHKERRQ(ierr); 1529 *n_global_subset = (PetscInt)PetscRealPart(globalsum); 1530 /* Fill global_vec with cumulative function for global numbering */ 1531 ierr = VecGetArray(global_vec,&array);CHKERRQ(ierr); 1532 ierr = VecGetLocalSize(global_vec,&s);CHKERRQ(ierr); 1533 nlocals = 0; 1534 first_index = -1; 1535 first_found = PETSC_FALSE; 1536 for (i=0;i<s;i++) { 1537 if (!first_found && PetscRealPart(array[i]) > 0.0) { 1538 first_found = PETSC_TRUE; 1539 first_index = i; 1540 } 1541 nlocals += (PetscInt)PetscRealPart(array[i]); 1542 } 1543 ierr = MPI_Gather(&nlocals,1,MPIU_INT,dof_sizes,1,MPIU_INT,0,comm);CHKERRQ(ierr); 1544 if (!rank_prec_comm) { 1545 dof_displs[0]=0; 1546 for (i=1;i<size_prec_comm;i++) { 1547 dof_displs[i] = dof_displs[i-1]+dof_sizes[i-1]; 1548 } 1549 } 1550 ierr = MPI_Scatter(dof_displs,1,MPIU_INT,&nlocals,1,MPIU_INT,0,comm);CHKERRQ(ierr); 1551 if (first_found) { 1552 array[first_index] += (PetscScalar)nlocals; 1553 old_index = first_index; 1554 for (i=first_index+1;i<s;i++) { 1555 if (PetscRealPart(array[i]) > 0.0) { 1556 array[i] += array[old_index]; 1557 old_index = i; 1558 } 1559 } 1560 } 1561 ierr = VecRestoreArray(global_vec,&array);CHKERRQ(ierr); 1562 ierr = VecSet(local_vec,0.0);CHKERRQ(ierr); 1563 ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1564 ierr = VecScatterEnd (scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1565 /* get global ordering of local dofs */ 1566 ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 1567 if (local_dofs_mult) { 1568 for (i=0;i<n_local_dofs;i++) { 1569 temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-local_dofs_mult[i]; 1570 } 1571 } else { 1572 for (i=0;i<n_local_dofs;i++) { 1573 temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-1; 1574 } 1575 } 1576 ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 1577 /* free workspace */ 1578 ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr); 1579 ierr = VecDestroy(&local_vec);CHKERRQ(ierr); 1580 ierr = VecDestroy(&global_vec);CHKERRQ(ierr); 1581 ierr = PetscFree(dof_sizes);CHKERRQ(ierr); 1582 ierr = PetscFree(dof_displs);CHKERRQ(ierr); 1583 /* return pointer to global ordering of local dofs */ 1584 *global_numbering_subset = temp_global_dofs; 1585 PetscFunctionReturn(0); 1586 } 1587 1588 #undef __FUNCT__ 1589 #define __FUNCT__ "PCBDDCOrthonormalizeVecs" 1590 PetscErrorCode PCBDDCOrthonormalizeVecs(PetscInt n, Vec vecs[]) 1591 { 1592 PetscInt i,j; 1593 PetscScalar *alphas; 1594 PetscErrorCode ierr; 1595 1596 PetscFunctionBegin; 1597 /* this implements stabilized Gram-Schmidt */ 1598 ierr = PetscMalloc(n*sizeof(PetscScalar),&alphas);CHKERRQ(ierr); 1599 for (i=0;i<n;i++) { 1600 ierr = VecNormalize(vecs[i],NULL);CHKERRQ(ierr); 1601 if (i<n) { ierr = VecMDot(vecs[i],n-i-1,&vecs[i+1],&alphas[i+1]);CHKERRQ(ierr); } 1602 for (j=i+1;j<n;j++) { ierr = VecAXPY(vecs[j],PetscConj(-alphas[j]),vecs[i]);CHKERRQ(ierr); } 1603 } 1604 ierr = PetscFree(alphas);CHKERRQ(ierr); 1605 PetscFunctionReturn(0); 1606 } 1607 1608 1609