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 = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 51 ierr = VecDestroy(&pcbddc->coarse_vec);CHKERRQ(ierr); 52 ierr = VecDestroy(&pcbddc->coarse_rhs);CHKERRQ(ierr); 53 ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr); 54 ierr = MatDestroy(&pcbddc->coarse_mat);CHKERRQ(ierr); 55 ierr = MatDestroy(&pcbddc->coarse_phi_B);CHKERRQ(ierr); 56 ierr = MatDestroy(&pcbddc->coarse_phi_D);CHKERRQ(ierr); 57 ierr = MatDestroy(&pcbddc->coarse_psi_B);CHKERRQ(ierr); 58 ierr = MatDestroy(&pcbddc->coarse_psi_D);CHKERRQ(ierr); 59 ierr = VecDestroy(&pcbddc->vec1_P);CHKERRQ(ierr); 60 ierr = VecDestroy(&pcbddc->vec1_C);CHKERRQ(ierr); 61 ierr = MatDestroy(&pcbddc->local_auxmat1);CHKERRQ(ierr); 62 ierr = MatDestroy(&pcbddc->local_auxmat2);CHKERRQ(ierr); 63 ierr = VecDestroy(&pcbddc->vec1_R);CHKERRQ(ierr); 64 ierr = VecDestroy(&pcbddc->vec2_R);CHKERRQ(ierr); 65 ierr = VecDestroy(&pcbddc->vec4_D);CHKERRQ(ierr); 66 ierr = VecScatterDestroy(&pcbddc->R_to_B);CHKERRQ(ierr); 67 ierr = VecScatterDestroy(&pcbddc->R_to_D);CHKERRQ(ierr); 68 ierr = VecScatterDestroy(&pcbddc->coarse_loc_to_glob);CHKERRQ(ierr); 69 ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr); 70 ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr); 71 ierr = PetscFree(pcbddc->local_primal_indices);CHKERRQ(ierr); 72 ierr = PetscFree(pcbddc->replicated_local_primal_indices);CHKERRQ(ierr); 73 ierr = PetscFree(pcbddc->replicated_local_primal_values);CHKERRQ(ierr); 74 ierr = PetscFree(pcbddc->local_primal_displacements);CHKERRQ(ierr); 75 ierr = PetscFree(pcbddc->local_primal_sizes);CHKERRQ(ierr); 76 PetscFunctionReturn(0); 77 } 78 79 #undef __FUNCT__ 80 #define __FUNCT__ "PCBDDCSolveSaddlePoint" 81 static PetscErrorCode PCBDDCSolveSaddlePoint(PC pc) 82 { 83 PetscErrorCode ierr; 84 PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 85 86 PetscFunctionBegin; 87 ierr = KSPSolve(pcbddc->ksp_R,pcbddc->vec1_R,pcbddc->vec2_R);CHKERRQ(ierr); 88 if (pcbddc->local_auxmat1) { 89 ierr = MatMult(pcbddc->local_auxmat1,pcbddc->vec2_R,pcbddc->vec1_C);CHKERRQ(ierr); 90 ierr = MatMultAdd(pcbddc->local_auxmat2,pcbddc->vec1_C,pcbddc->vec2_R,pcbddc->vec2_R);CHKERRQ(ierr); 91 } 92 PetscFunctionReturn(0); 93 } 94 95 #undef __FUNCT__ 96 #define __FUNCT__ "PCBDDCApplyInterfacePreconditioner" 97 PetscErrorCode PCBDDCApplyInterfacePreconditioner(PC pc) 98 { 99 PetscErrorCode ierr; 100 PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 101 PC_IS* pcis = (PC_IS*) (pc->data); 102 const PetscScalar zero = 0.0; 103 104 PetscFunctionBegin; 105 /* Application of PHI^T (or PSI^T) */ 106 if (pcbddc->coarse_psi_B) { 107 ierr = MatMultTranspose(pcbddc->coarse_psi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr); 108 if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_psi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); } 109 } else { 110 ierr = MatMultTranspose(pcbddc->coarse_phi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr); 111 if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_phi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); } 112 } 113 /* Scatter data of coarse_rhs */ 114 if (pcbddc->coarse_rhs) { ierr = VecSet(pcbddc->coarse_rhs,zero);CHKERRQ(ierr); } 115 ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 116 117 /* Local solution on R nodes */ 118 ierr = VecSet(pcbddc->vec1_R,zero);CHKERRQ(ierr); 119 ierr = VecScatterBegin(pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 120 ierr = VecScatterEnd (pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 121 if (pcbddc->inexact_prec_type) { 122 ierr = VecScatterBegin(pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 123 ierr = VecScatterEnd (pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 124 } 125 ierr = PCBDDCSolveSaddlePoint(pc);CHKERRQ(ierr); 126 ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr); 127 ierr = VecScatterBegin(pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 128 ierr = VecScatterEnd (pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 129 if (pcbddc->inexact_prec_type) { 130 ierr = VecScatterBegin(pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 131 ierr = VecScatterEnd (pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 132 } 133 134 /* Coarse solution */ 135 ierr = PCBDDCScatterCoarseDataEnd(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 136 if (pcbddc->coarse_rhs) { /* TODO remove null space when doing multilevel */ 137 ierr = KSPSolve(pcbddc->coarse_ksp,pcbddc->coarse_rhs,pcbddc->coarse_vec);CHKERRQ(ierr); 138 } 139 ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 140 ierr = PCBDDCScatterCoarseDataEnd (pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 141 142 /* Sum contributions from two levels */ 143 ierr = MatMultAdd(pcbddc->coarse_phi_B,pcbddc->vec1_P,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr); 144 if (pcbddc->inexact_prec_type) { ierr = MatMultAdd(pcbddc->coarse_phi_D,pcbddc->vec1_P,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 145 PetscFunctionReturn(0); 146 } 147 148 #undef __FUNCT__ 149 #define __FUNCT__ "PCBDDCScatterCoarseDataBegin" 150 PetscErrorCode PCBDDCScatterCoarseDataBegin(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode) 151 { 152 PetscErrorCode ierr; 153 PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 154 155 PetscFunctionBegin; 156 switch (pcbddc->coarse_communications_type) { 157 case SCATTERS_BDDC: 158 ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr); 159 break; 160 case GATHERS_BDDC: 161 break; 162 } 163 PetscFunctionReturn(0); 164 } 165 166 #undef __FUNCT__ 167 #define __FUNCT__ "PCBDDCScatterCoarseDataEnd" 168 PetscErrorCode PCBDDCScatterCoarseDataEnd(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode) 169 { 170 PetscErrorCode ierr; 171 PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 172 PetscScalar* array_to; 173 PetscScalar* array_from; 174 MPI_Comm comm; 175 PetscInt i; 176 177 PetscFunctionBegin; 178 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 179 switch (pcbddc->coarse_communications_type) { 180 case SCATTERS_BDDC: 181 ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr); 182 break; 183 case GATHERS_BDDC: 184 if (vec_from) { 185 ierr = VecGetArray(vec_from,&array_from);CHKERRQ(ierr); 186 } 187 if (vec_to) { 188 ierr = VecGetArray(vec_to,&array_to);CHKERRQ(ierr); 189 } 190 switch(pcbddc->coarse_problem_type){ 191 case SEQUENTIAL_BDDC: 192 if (smode == SCATTER_FORWARD) { 193 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); 194 if (vec_to) { 195 if (imode == ADD_VALUES) { 196 for (i=0;i<pcbddc->replicated_primal_size;i++) { 197 array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i]; 198 } 199 } else { 200 for (i=0;i<pcbddc->replicated_primal_size;i++) { 201 array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i]; 202 } 203 } 204 } 205 } else { 206 if (vec_from) { 207 if (imode == ADD_VALUES) { 208 MPI_Comm vec_from_comm; 209 ierr = PetscObjectGetComm((PetscObject)(vec_from),&vec_from_comm);CHKERRQ(ierr); 210 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); 211 } 212 for (i=0;i<pcbddc->replicated_primal_size;i++) { 213 pcbddc->replicated_local_primal_values[i]=array_from[pcbddc->replicated_local_primal_indices[i]]; 214 } 215 } 216 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); 217 } 218 break; 219 case REPLICATED_BDDC: 220 if (smode == SCATTER_FORWARD) { 221 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); 222 if (imode == ADD_VALUES) { 223 for (i=0;i<pcbddc->replicated_primal_size;i++) { 224 array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i]; 225 } 226 } else { 227 for (i=0;i<pcbddc->replicated_primal_size;i++) { 228 array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i]; 229 } 230 } 231 } else { /* no communications needed for SCATTER_REVERSE since needed data is already present */ 232 if (imode == ADD_VALUES) { 233 for (i=0;i<pcbddc->local_primal_size;i++) { 234 array_to[i]+=array_from[pcbddc->local_primal_indices[i]]; 235 } 236 } else { 237 for (i=0;i<pcbddc->local_primal_size;i++) { 238 array_to[i]=array_from[pcbddc->local_primal_indices[i]]; 239 } 240 } 241 } 242 break; 243 case MULTILEVEL_BDDC: 244 break; 245 case PARALLEL_BDDC: 246 break; 247 } 248 if (vec_from) { 249 ierr = VecRestoreArray(vec_from,&array_from);CHKERRQ(ierr); 250 } 251 if (vec_to) { 252 ierr = VecRestoreArray(vec_to,&array_to);CHKERRQ(ierr); 253 } 254 break; 255 } 256 PetscFunctionReturn(0); 257 } 258 259 #undef __FUNCT__ 260 #define __FUNCT__ "PCBDDCConstraintsSetUp" 261 PetscErrorCode PCBDDCConstraintsSetUp(PC pc) 262 { 263 PetscErrorCode ierr; 264 PC_IS* pcis = (PC_IS*)(pc->data); 265 PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 266 Mat_IS *matis = (Mat_IS*)pc->pmat->data; 267 PetscInt *nnz,*is_indices; 268 PetscScalar *temp_quadrature_constraint; 269 PetscInt *temp_indices,*temp_indices_to_constraint,*temp_indices_to_constraint_B,*local_to_B; 270 PetscInt local_primal_size,i,j,k,total_counts,max_size_of_constraint; 271 PetscInt n_vertices,size_of_constraint; 272 PetscScalar quad_value; 273 PetscBool nnsp_has_cnst=PETSC_FALSE,use_nnsp_true=pcbddc->use_nnsp_true; 274 PetscInt nnsp_size=0,nnsp_addone=0,temp_constraints,temp_start_ptr,n_ISForFaces,n_ISForEdges; 275 IS *used_IS,ISForVertices,*ISForFaces,*ISForEdges; 276 MatType impMatType=MATSEQAIJ; 277 PetscBLASInt Bs,Bt,lwork,lierr; 278 PetscReal tol=1.0e-8; 279 MatNullSpace nearnullsp; 280 const Vec *nearnullvecs; 281 Vec *localnearnullsp; 282 PetscScalar *work,*temp_basis,*array_vector,*correlation_mat; 283 PetscReal *rwork,*singular_vals; 284 PetscBLASInt Bone=1,*ipiv; 285 Vec temp_vec; 286 Mat temp_mat; 287 KSP temp_ksp; 288 PC temp_pc; 289 PetscInt s,start_constraint,dual_dofs; 290 PetscBool compute_submatrix,useksp=PETSC_FALSE; 291 PetscInt *aux_primal_permutation,*aux_primal_numbering; 292 PetscBool boolforchange,*change_basis; 293 /* some ugly conditional declarations */ 294 #if defined(PETSC_MISSING_LAPACK_GESVD) 295 PetscScalar one=1.0,zero=0.0; 296 PetscInt ii; 297 PetscScalar *singular_vectors; 298 PetscBLASInt *iwork,*ifail; 299 PetscReal dummy_real,abs_tol; 300 PetscBLASInt eigs_found; 301 #endif 302 PetscBLASInt dummy_int; 303 PetscScalar dummy_scalar; 304 PetscBool used_vertex,get_faces,get_edges,get_vertices; 305 306 PetscFunctionBegin; 307 /* Get index sets for faces, edges and vertices from graph */ 308 get_faces = PETSC_TRUE; 309 get_edges = PETSC_TRUE; 310 get_vertices = PETSC_TRUE; 311 if (pcbddc->vertices_flag) { 312 get_faces = PETSC_FALSE; 313 get_edges = PETSC_FALSE; 314 } 315 if (pcbddc->constraints_flag) { 316 get_vertices = PETSC_FALSE; 317 } 318 if (pcbddc->faces_flag) { 319 get_edges = PETSC_FALSE; 320 } 321 if (pcbddc->edges_flag) { 322 get_faces = PETSC_FALSE; 323 } 324 /* default */ 325 if (!get_faces && !get_edges && !get_vertices) { 326 get_vertices = PETSC_TRUE; 327 } 328 ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,get_faces,get_edges,get_vertices,&n_ISForFaces,&ISForFaces,&n_ISForEdges,&ISForEdges,&ISForVertices); 329 if (pcbddc->dbg_flag) { 330 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr); 331 i = 0; 332 if (ISForVertices) { 333 ierr = ISGetSize(ISForVertices,&i);CHKERRQ(ierr); 334 } 335 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate vertices\n",PetscGlobalRank,i);CHKERRQ(ierr); 336 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate edges\n",PetscGlobalRank,n_ISForEdges);CHKERRQ(ierr); 337 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate faces\n",PetscGlobalRank,n_ISForFaces);CHKERRQ(ierr); 338 ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 339 } 340 /* check if near null space is attached to global mat */ 341 ierr = MatGetNearNullSpace(pc->pmat,&nearnullsp);CHKERRQ(ierr); 342 if (nearnullsp) { 343 ierr = MatNullSpaceGetVecs(nearnullsp,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr); 344 } else { /* if near null space is not provided it uses constants */ 345 nnsp_has_cnst = PETSC_TRUE; 346 use_nnsp_true = PETSC_TRUE; 347 } 348 if (nnsp_has_cnst) { 349 nnsp_addone = 1; 350 } 351 /* 352 Evaluate maximum storage size needed by the procedure 353 - temp_indices will contain start index of each constraint stored as follows 354 - temp_indices_to_constraint [temp_indices[i],...,temp[indices[i+1]-1] will contain the indices (in local numbering) on which the constraint acts 355 - 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 356 - temp_quadrature_constraint [temp_indices[i],...,temp[indices[i+1]-1] will contain the scalars representing the constraint itself 357 */ 358 total_counts = n_ISForFaces+n_ISForEdges; 359 total_counts *= (nnsp_addone+nnsp_size); 360 n_vertices = 0; 361 if (ISForVertices) { 362 ierr = ISGetSize(ISForVertices,&n_vertices);CHKERRQ(ierr); 363 } 364 total_counts += n_vertices; 365 ierr = PetscMalloc((total_counts+1)*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr); 366 ierr = PetscMalloc((total_counts+1)*sizeof(PetscBool),&change_basis);CHKERRQ(ierr); 367 total_counts = 0; 368 max_size_of_constraint = 0; 369 for (i=0;i<n_ISForEdges+n_ISForFaces;i++) { 370 if (i<n_ISForEdges) { 371 used_IS = &ISForEdges[i]; 372 } else { 373 used_IS = &ISForFaces[i-n_ISForEdges]; 374 } 375 ierr = ISGetSize(*used_IS,&j);CHKERRQ(ierr); 376 total_counts += j; 377 max_size_of_constraint = PetscMax(j,max_size_of_constraint); 378 } 379 total_counts *= (nnsp_addone+nnsp_size); 380 total_counts += n_vertices; 381 ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&temp_quadrature_constraint);CHKERRQ(ierr); 382 ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint);CHKERRQ(ierr); 383 ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint_B);CHKERRQ(ierr); 384 ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&local_to_B);CHKERRQ(ierr); 385 ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 386 for (i=0;i<pcis->n;i++) { 387 local_to_B[i]=-1; 388 } 389 for (i=0;i<pcis->n_B;i++) { 390 local_to_B[is_indices[i]]=i; 391 } 392 ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 393 394 /* First we issue queries to allocate optimal workspace for LAPACKgesvd or LAPACKsyev/LAPACKheev */ 395 rwork = 0; 396 work = 0; 397 singular_vals = 0; 398 temp_basis = 0; 399 correlation_mat = 0; 400 if (!pcbddc->use_nnsp_true) { 401 PetscScalar temp_work; 402 #if defined(PETSC_MISSING_LAPACK_GESVD) 403 /* POD */ 404 PetscInt max_n; 405 max_n = nnsp_addone+nnsp_size; 406 /* using some techniques borrowed from Proper Orthogonal Decomposition */ 407 ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&correlation_mat);CHKERRQ(ierr); 408 ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&singular_vectors);CHKERRQ(ierr); 409 ierr = PetscMalloc(max_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr); 410 ierr = PetscMalloc(max_size_of_constraint*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr); 411 #if defined(PETSC_USE_COMPLEX) 412 ierr = PetscMalloc(3*max_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr); 413 #endif 414 ierr = PetscMalloc(5*max_n*sizeof(PetscBLASInt),&iwork);CHKERRQ(ierr); 415 ierr = PetscMalloc(max_n*sizeof(PetscBLASInt),&ifail);CHKERRQ(ierr); 416 /* now we evaluate the optimal workspace using query with lwork=-1 */ 417 ierr = PetscBLASIntCast(max_n,&Bt);CHKERRQ(ierr); 418 lwork=-1; 419 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 420 #if !defined(PETSC_USE_COMPLEX) 421 abs_tol=1.e-8; 422 /* LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,&lierr); */ 423 LAPACKsyevx_("V","A","U",&Bt,correlation_mat,&Bt,&dummy_real,&dummy_real,&dummy_int,&dummy_int, 424 &abs_tol,&eigs_found,singular_vals,singular_vectors,&Bt,&temp_work,&lwork,iwork,ifail,&lierr); 425 #else 426 /* LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,rwork,&lierr); */ 427 /* LAPACK call is missing here! TODO */ 428 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1"); 429 #endif 430 if ( lierr ) { 431 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYEVX Lapack routine %d",(int)lierr); 432 } 433 ierr = PetscFPTrapPop();CHKERRQ(ierr); 434 #else /* on missing GESVD */ 435 /* SVD */ 436 PetscInt max_n,min_n; 437 max_n = max_size_of_constraint; 438 min_n = nnsp_addone+nnsp_size; 439 if (max_size_of_constraint < ( nnsp_addone+nnsp_size ) ) { 440 min_n = max_size_of_constraint; 441 max_n = nnsp_addone+nnsp_size; 442 } 443 ierr = PetscMalloc(min_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr); 444 #if defined(PETSC_USE_COMPLEX) 445 ierr = PetscMalloc(5*min_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr); 446 #endif 447 /* now we evaluate the optimal workspace using query with lwork=-1 */ 448 lwork=-1; 449 ierr = PetscBLASIntCast(max_n,&Bs);CHKERRQ(ierr); 450 ierr = PetscBLASIntCast(min_n,&Bt);CHKERRQ(ierr); 451 dummy_int = Bs; 452 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 453 #if !defined(PETSC_USE_COMPLEX) 454 LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Bs,singular_vals, 455 &dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,&temp_work,&lwork,&lierr); 456 #else 457 LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Bs,singular_vals, 458 &dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,&temp_work,&lwork,rwork,&lierr); 459 #endif 460 if ( lierr ) { 461 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SVD Lapack routine %d",(int)lierr); 462 } 463 ierr = PetscFPTrapPop();CHKERRQ(ierr); 464 #endif 465 /* Allocate optimal workspace */ 466 ierr = PetscBLASIntCast((PetscInt)PetscRealPart(temp_work),&lwork);CHKERRQ(ierr); 467 total_counts = (PetscInt)lwork; 468 ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&work);CHKERRQ(ierr); 469 } 470 /* get local part of global near null space vectors */ 471 ierr = PetscMalloc(nnsp_size*sizeof(Vec),&localnearnullsp);CHKERRQ(ierr); 472 for (k=0;k<nnsp_size;k++) { 473 ierr = VecDuplicate(pcis->vec1_N,&localnearnullsp[k]);CHKERRQ(ierr); 474 ierr = VecScatterBegin(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 475 ierr = VecScatterEnd(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 476 } 477 /* Now we can loop on constraining sets */ 478 total_counts = 0; 479 temp_indices[0] = 0; 480 /* vertices */ 481 if (ISForVertices) { 482 ierr = ISGetIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); 483 if (nnsp_has_cnst) { /* consider all vertices */ 484 for (i=0;i<n_vertices;i++) { 485 temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i]; 486 temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]]; 487 temp_quadrature_constraint[temp_indices[total_counts]]=1.0; 488 temp_indices[total_counts+1]=temp_indices[total_counts]+1; 489 change_basis[total_counts]=PETSC_FALSE; 490 total_counts++; 491 } 492 } else { /* consider vertices for which exist at least a localnearnullsp which is not null there */ 493 for (i=0;i<n_vertices;i++) { 494 used_vertex=PETSC_FALSE; 495 k=0; 496 while (!used_vertex && k<nnsp_size) { 497 ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 498 if (PetscAbsScalar(array_vector[is_indices[i]])>0.0) { 499 temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i]; 500 temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]]; 501 temp_quadrature_constraint[temp_indices[total_counts]]=1.0; 502 temp_indices[total_counts+1]=temp_indices[total_counts]+1; 503 change_basis[total_counts]=PETSC_FALSE; 504 total_counts++; 505 used_vertex=PETSC_TRUE; 506 } 507 ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 508 k++; 509 } 510 } 511 } 512 ierr = ISRestoreIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); 513 n_vertices = total_counts; 514 } 515 /* edges and faces */ 516 for (i=0;i<n_ISForEdges+n_ISForFaces;i++) { 517 if (i<n_ISForEdges) { 518 used_IS = &ISForEdges[i]; 519 boolforchange = pcbddc->use_change_of_basis; 520 } else { 521 used_IS = &ISForFaces[i-n_ISForEdges]; 522 boolforchange = pcbddc->use_change_on_faces; 523 } 524 temp_constraints = 0; /* zero the number of constraints I have on this conn comp */ 525 temp_start_ptr = total_counts; /* need to know the starting index of constraints stored */ 526 ierr = ISGetSize(*used_IS,&size_of_constraint);CHKERRQ(ierr); 527 ierr = ISGetIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 528 /* HACK: change of basis should not performed on local periodic nodes */ 529 if (pcbddc->mat_graph->mirrors && pcbddc->mat_graph->mirrors[is_indices[0]]) { 530 boolforchange = PETSC_FALSE; 531 } 532 if (nnsp_has_cnst) { 533 temp_constraints++; 534 quad_value = (PetscScalar)(1.0/PetscSqrtReal((PetscReal)size_of_constraint)); 535 for (j=0;j<size_of_constraint;j++) { 536 temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j]; 537 temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]]; 538 temp_quadrature_constraint[temp_indices[total_counts]+j]=quad_value; 539 } 540 temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint; /* store new starting point */ 541 change_basis[total_counts]=boolforchange; 542 total_counts++; 543 } 544 for (k=0;k<nnsp_size;k++) { 545 ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 546 for (j=0;j<size_of_constraint;j++) { 547 temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j]; 548 temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]]; 549 temp_quadrature_constraint[temp_indices[total_counts]+j]=array_vector[is_indices[j]]; 550 } 551 ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 552 quad_value = 1.0; 553 if (use_nnsp_true) { /* check if array is null on the connected component in case use_nnsp_true has been requested */ 554 ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr); 555 quad_value = BLASasum_(&Bs,&temp_quadrature_constraint[temp_indices[total_counts]],&Bone); 556 } 557 if (quad_value > 0.0) { /* keep indices and values */ 558 temp_constraints++; 559 temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint; /* store new starting point */ 560 change_basis[total_counts]=boolforchange; 561 total_counts++; 562 } 563 } 564 ierr = ISRestoreIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 565 /* perform SVD on the constraint if use_nnsp_true has not be requested by the user */ 566 if (!use_nnsp_true) { 567 ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr); 568 ierr = PetscBLASIntCast(temp_constraints,&Bt);CHKERRQ(ierr); 569 570 #if defined(PETSC_MISSING_LAPACK_GESVD) 571 ierr = PetscMemzero(correlation_mat,Bt*Bt*sizeof(PetscScalar));CHKERRQ(ierr); 572 /* Store upper triangular part of correlation matrix */ 573 for (j=0;j<temp_constraints;j++) { 574 for (k=0;k<j+1;k++) { 575 correlation_mat[j*temp_constraints+k]=BLASdot_(&Bs,&temp_quadrature_constraint[temp_indices[temp_start_ptr+j]],&Bone, 576 &temp_quadrature_constraint[temp_indices[temp_start_ptr+k]],&Bone); 577 578 } 579 } 580 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 581 #if !defined(PETSC_USE_COMPLEX) 582 /* LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,work,&lwork,&lierr); */ 583 LAPACKsyevx_("V","A","U",&Bt,correlation_mat,&Bt,&dummy_real,&dummy_real,&dummy_int,&dummy_int, 584 &abs_tol,&eigs_found,singular_vals,singular_vectors,&Bt,work,&lwork,iwork,ifail,&lierr); 585 #else 586 /* LAPACK call is missing here! TODO */ 587 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1"); 588 #endif 589 if (lierr) { 590 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYEVX Lapack routine %d",(int)lierr); 591 } 592 ierr = PetscFPTrapPop();CHKERRQ(ierr); 593 /* retain eigenvalues greater than tol: note that lapack SYEV gives eigs in ascending order */ 594 j=0; 595 while (j < Bt && singular_vals[j] < tol) j++; 596 total_counts=total_counts-j; 597 if (j<temp_constraints) { 598 for (k=j;k<Bt;k++) { 599 singular_vals[k]=1.0/PetscSqrtReal(singular_vals[k]); 600 } 601 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 602 BLASgemm_("N","N",&Bs,&Bt,&Bt,&one,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,correlation_mat,&Bt,&zero,temp_basis,&Bs); 603 ierr = PetscFPTrapPop();CHKERRQ(ierr); 604 /* copy POD basis into used quadrature memory */ 605 for (k=0;k<Bt-j;k++) { 606 for (ii=0;ii<size_of_constraint;ii++) { 607 temp_quadrature_constraint[temp_indices[temp_start_ptr+k]+ii]=singular_vals[Bt-1-k]*temp_basis[(Bt-1-k)*size_of_constraint+ii]; 608 } 609 } 610 } 611 612 #else /* on missing GESVD */ 613 PetscInt min_n = temp_constraints; 614 if (min_n > size_of_constraint) min_n = size_of_constraint; 615 dummy_int = Bs; 616 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 617 #if !defined(PETSC_USE_COMPLEX) 618 LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,singular_vals, 619 &dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,work,&lwork,&lierr); 620 #else 621 LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,singular_vals, 622 &dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,work,&lwork,rwork,&lierr); 623 #endif 624 if (lierr) { 625 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SVD Lapack routine %d",(int)lierr); 626 } 627 ierr = PetscFPTrapPop();CHKERRQ(ierr); 628 /* retain eigenvalues greater than tol: note that lapack SVD gives eigs in descending order */ 629 j = 0; 630 while (j < min_n && singular_vals[min_n-j-1] < tol) j++; 631 total_counts = total_counts-(PetscInt)Bt+(min_n-j); 632 #endif 633 } 634 } 635 /* free index sets of faces, edges and vertices */ 636 for (i=0;i<n_ISForFaces;i++) { 637 ierr = ISDestroy(&ISForFaces[i]);CHKERRQ(ierr); 638 } 639 ierr = PetscFree(ISForFaces);CHKERRQ(ierr); 640 for (i=0;i<n_ISForEdges;i++) { 641 ierr = ISDestroy(&ISForEdges[i]);CHKERRQ(ierr); 642 } 643 ierr = PetscFree(ISForEdges);CHKERRQ(ierr); 644 ierr = ISDestroy(&ISForVertices);CHKERRQ(ierr); 645 646 /* set quantities in pcbddc data structure */ 647 /* n_vertices defines the number of point primal dofs */ 648 /* n_constraints defines the number of averages (they can be point primal dofs if change of basis is requested) */ 649 local_primal_size = total_counts; 650 pcbddc->n_vertices = n_vertices; 651 pcbddc->n_constraints = total_counts-n_vertices; 652 pcbddc->local_primal_size = local_primal_size; 653 654 /* Create constraint matrix */ 655 /* The constraint matrix is used to compute the l2g map of primal dofs */ 656 /* so we need to set it up properly either with or without change of basis */ 657 ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ConstraintMatrix);CHKERRQ(ierr); 658 ierr = MatSetType(pcbddc->ConstraintMatrix,impMatType);CHKERRQ(ierr); 659 ierr = MatSetSizes(pcbddc->ConstraintMatrix,local_primal_size,pcis->n,local_primal_size,pcis->n);CHKERRQ(ierr); 660 /* compute a local numbering of constraints : vertices first then constraints */ 661 ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr); 662 ierr = VecGetArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr); 663 ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_numbering);CHKERRQ(ierr); 664 ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_permutation);CHKERRQ(ierr); 665 total_counts=0; 666 /* find vertices: subdomain corners plus dofs with basis changed */ 667 for (i=0;i<local_primal_size;i++) { 668 size_of_constraint=temp_indices[i+1]-temp_indices[i]; 669 if (change_basis[i] || size_of_constraint == 1) { 670 k=0; 671 while(k < size_of_constraint && array_vector[temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1]] != 0.0) { 672 k=k+1; 673 } 674 j=temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1]; 675 array_vector[j] = 1.0; 676 aux_primal_numbering[total_counts]=j; 677 aux_primal_permutation[total_counts]=total_counts; 678 total_counts++; 679 } 680 } 681 ierr = VecRestoreArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr); 682 /* permute indices in order to have a sorted set of vertices */ 683 ierr = PetscSortIntWithPermutation(total_counts,aux_primal_numbering,aux_primal_permutation); 684 /* nonzero structure */ 685 ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 686 for (i=0;i<total_counts;i++) { 687 nnz[i]=1; 688 } 689 j=total_counts; 690 for (i=n_vertices;i<local_primal_size;i++) { 691 if (!change_basis[i]) { 692 nnz[j]=temp_indices[i+1]-temp_indices[i]; 693 j++; 694 } 695 } 696 ierr = MatSeqAIJSetPreallocation(pcbddc->ConstraintMatrix,0,nnz);CHKERRQ(ierr); 697 ierr = PetscFree(nnz);CHKERRQ(ierr); 698 /* set values in constraint matrix */ 699 for (i=0;i<total_counts;i++) { 700 j = aux_primal_permutation[i]; 701 k = aux_primal_numbering[j]; 702 ierr = MatSetValue(pcbddc->ConstraintMatrix,i,k,1.0,INSERT_VALUES);CHKERRQ(ierr); 703 } 704 for (i=n_vertices;i<local_primal_size;i++) { 705 if (!change_basis[i]) { 706 size_of_constraint=temp_indices[i+1]-temp_indices[i]; 707 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); 708 total_counts++; 709 } 710 } 711 ierr = PetscFree(aux_primal_numbering);CHKERRQ(ierr); 712 ierr = PetscFree(aux_primal_permutation);CHKERRQ(ierr); 713 /* assembling */ 714 ierr = MatAssemblyBegin(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 715 ierr = MatAssemblyEnd(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 716 717 /* Create matrix for change of basis. We don't need it in case pcbddc->use_change_of_basis is FALSE */ 718 if (pcbddc->use_change_of_basis) { 719 ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr); 720 ierr = MatSetType(pcbddc->ChangeOfBasisMatrix,impMatType);CHKERRQ(ierr); 721 ierr = MatSetSizes(pcbddc->ChangeOfBasisMatrix,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);CHKERRQ(ierr); 722 /* work arrays */ 723 /* we need to reuse these arrays, so we free them */ 724 ierr = PetscFree(temp_basis);CHKERRQ(ierr); 725 ierr = PetscFree(work);CHKERRQ(ierr); 726 ierr = PetscMalloc(pcis->n_B*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 727 ierr = PetscMalloc((nnsp_addone+nnsp_size)*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr); 728 ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscScalar),&work);CHKERRQ(ierr); 729 ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscBLASInt),&ipiv);CHKERRQ(ierr); 730 for (i=0;i<pcis->n_B;i++) { 731 nnz[i]=1; 732 } 733 /* Overestimated nonzeros per row */ 734 k=1; 735 for (i=pcbddc->n_vertices;i<local_primal_size;i++) { 736 if (change_basis[i]) { 737 size_of_constraint = temp_indices[i+1]-temp_indices[i]; 738 if (k < size_of_constraint) { 739 k = size_of_constraint; 740 } 741 for (j=0;j<size_of_constraint;j++) { 742 nnz[temp_indices_to_constraint_B[temp_indices[i]+j]] = size_of_constraint; 743 } 744 } 745 } 746 ierr = MatSeqAIJSetPreallocation(pcbddc->ChangeOfBasisMatrix,0,nnz);CHKERRQ(ierr); 747 ierr = PetscFree(nnz);CHKERRQ(ierr); 748 /* Temporary array to store indices */ 749 ierr = PetscMalloc(k*sizeof(PetscInt),&is_indices);CHKERRQ(ierr); 750 /* Set initial identity in the matrix */ 751 for (i=0;i<pcis->n_B;i++) { 752 ierr = MatSetValue(pcbddc->ChangeOfBasisMatrix,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr); 753 } 754 /* Now we loop on the constraints which need a change of basis */ 755 /* Change of basis matrix is evaluated as the FIRST APPROACH in */ 756 /* Klawonn and Widlund, Dual-primal FETI-DP methods for linear elasticity, (6.2.1) */ 757 temp_constraints = 0; 758 if (pcbddc->n_vertices < local_primal_size) { 759 temp_start_ptr = temp_indices_to_constraint_B[temp_indices[pcbddc->n_vertices]]; 760 } 761 for (i=pcbddc->n_vertices;i<local_primal_size;i++) { 762 if (change_basis[i]) { 763 compute_submatrix = PETSC_FALSE; 764 useksp = PETSC_FALSE; 765 if (temp_start_ptr == temp_indices_to_constraint_B[temp_indices[i]]) { 766 temp_constraints++; 767 if (i == local_primal_size -1 || temp_start_ptr != temp_indices_to_constraint_B[temp_indices[i+1]]) { 768 compute_submatrix = PETSC_TRUE; 769 } 770 } 771 if (compute_submatrix) { 772 if (temp_constraints > 1 || pcbddc->use_nnsp_true) { 773 useksp = PETSC_TRUE; 774 } 775 size_of_constraint = temp_indices[i+1]-temp_indices[i]; 776 if (useksp) { /* experimental TODO: reuse KSP and MAT instead of creating them each time */ 777 ierr = MatCreate(PETSC_COMM_SELF,&temp_mat);CHKERRQ(ierr); 778 ierr = MatSetType(temp_mat,impMatType);CHKERRQ(ierr); 779 ierr = MatSetSizes(temp_mat,size_of_constraint,size_of_constraint,size_of_constraint,size_of_constraint);CHKERRQ(ierr); 780 ierr = MatSeqAIJSetPreallocation(temp_mat,size_of_constraint,NULL);CHKERRQ(ierr); 781 } 782 /* First _size_of_constraint-temp_constraints_ columns */ 783 dual_dofs = size_of_constraint-temp_constraints; 784 start_constraint = i+1-temp_constraints; 785 for (s=0;s<dual_dofs;s++) { 786 is_indices[0] = s; 787 for (j=0;j<temp_constraints;j++) { 788 for (k=0;k<temp_constraints;k++) { 789 temp_basis[j*temp_constraints+k]=temp_quadrature_constraint[temp_indices[start_constraint+k]+s+j+1]; 790 } 791 work[j]=-temp_quadrature_constraint[temp_indices[start_constraint+j]+s]; 792 is_indices[j+1]=s+j+1; 793 } 794 Bt = temp_constraints; 795 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 796 LAPACKgesv_(&Bt,&Bone,temp_basis,&Bt,ipiv,work,&Bt,&lierr); 797 if ( lierr ) { 798 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GESV Lapack routine %d",(int)lierr); 799 } 800 ierr = PetscFPTrapPop();CHKERRQ(ierr); 801 j = temp_indices_to_constraint_B[temp_indices[start_constraint]+s]; 802 ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,temp_constraints,&temp_indices_to_constraint_B[temp_indices[start_constraint]+s+1],1,&j,work,INSERT_VALUES);CHKERRQ(ierr); 803 if (useksp) { 804 /* temp mat with transposed rows and columns */ 805 ierr = MatSetValues(temp_mat,1,&s,temp_constraints,&is_indices[1],work,INSERT_VALUES);CHKERRQ(ierr); 806 ierr = MatSetValue(temp_mat,is_indices[0],is_indices[0],1.0,INSERT_VALUES);CHKERRQ(ierr); 807 } 808 } 809 if (useksp) { 810 /* last rows of temp_mat */ 811 for (j=0;j<size_of_constraint;j++) { 812 is_indices[j] = j; 813 } 814 for (s=0;s<temp_constraints;s++) { 815 k = s + dual_dofs; 816 ierr = MatSetValues(temp_mat,1,&k,size_of_constraint,is_indices,&temp_quadrature_constraint[temp_indices[start_constraint+s]],INSERT_VALUES);CHKERRQ(ierr); 817 } 818 ierr = MatAssemblyBegin(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 819 ierr = MatAssemblyEnd(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 820 ierr = MatGetVecs(temp_mat,&temp_vec,NULL);CHKERRQ(ierr); 821 ierr = KSPCreate(PETSC_COMM_SELF,&temp_ksp);CHKERRQ(ierr); 822 ierr = KSPSetOperators(temp_ksp,temp_mat,temp_mat,SAME_PRECONDITIONER);CHKERRQ(ierr); 823 ierr = KSPSetType(temp_ksp,KSPPREONLY);CHKERRQ(ierr); 824 ierr = KSPGetPC(temp_ksp,&temp_pc);CHKERRQ(ierr); 825 ierr = PCSetType(temp_pc,PCLU);CHKERRQ(ierr); 826 ierr = KSPSetUp(temp_ksp);CHKERRQ(ierr); 827 for (s=0;s<temp_constraints;s++) { 828 ierr = VecSet(temp_vec,0.0);CHKERRQ(ierr); 829 ierr = VecSetValue(temp_vec,s+dual_dofs,1.0,INSERT_VALUES);CHKERRQ(ierr); 830 ierr = VecAssemblyBegin(temp_vec);CHKERRQ(ierr); 831 ierr = VecAssemblyEnd(temp_vec);CHKERRQ(ierr); 832 ierr = KSPSolve(temp_ksp,temp_vec,temp_vec);CHKERRQ(ierr); 833 ierr = VecGetArray(temp_vec,&array_vector);CHKERRQ(ierr); 834 j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1]; 835 /* last columns of change of basis matrix associated to new primal dofs */ 836 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); 837 ierr = VecRestoreArray(temp_vec,&array_vector);CHKERRQ(ierr); 838 } 839 ierr = MatDestroy(&temp_mat);CHKERRQ(ierr); 840 ierr = KSPDestroy(&temp_ksp);CHKERRQ(ierr); 841 ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 842 } else { 843 /* last columns of change of basis matrix associated to new primal dofs */ 844 for (s=0;s<temp_constraints;s++) { 845 j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1]; 846 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); 847 } 848 } 849 /* prepare for the next cycle */ 850 temp_constraints = 0; 851 if (i != local_primal_size -1 ) { 852 temp_start_ptr = temp_indices_to_constraint_B[temp_indices[i+1]]; 853 } 854 } 855 } 856 } 857 /* assembling */ 858 ierr = MatAssemblyBegin(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 859 ierr = MatAssemblyEnd(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 860 ierr = PetscFree(ipiv);CHKERRQ(ierr); 861 ierr = PetscFree(is_indices);CHKERRQ(ierr); 862 } 863 /* free workspace no longer needed */ 864 ierr = PetscFree(rwork);CHKERRQ(ierr); 865 ierr = PetscFree(work);CHKERRQ(ierr); 866 ierr = PetscFree(temp_basis);CHKERRQ(ierr); 867 ierr = PetscFree(singular_vals);CHKERRQ(ierr); 868 ierr = PetscFree(correlation_mat);CHKERRQ(ierr); 869 ierr = PetscFree(temp_indices);CHKERRQ(ierr); 870 ierr = PetscFree(change_basis);CHKERRQ(ierr); 871 ierr = PetscFree(temp_indices_to_constraint);CHKERRQ(ierr); 872 ierr = PetscFree(temp_indices_to_constraint_B);CHKERRQ(ierr); 873 ierr = PetscFree(local_to_B);CHKERRQ(ierr); 874 ierr = PetscFree(temp_quadrature_constraint);CHKERRQ(ierr); 875 #if defined(PETSC_MISSING_LAPACK_GESVD) 876 ierr = PetscFree(iwork);CHKERRQ(ierr); 877 ierr = PetscFree(ifail);CHKERRQ(ierr); 878 ierr = PetscFree(singular_vectors);CHKERRQ(ierr); 879 #endif 880 for (k=0;k<nnsp_size;k++) { 881 ierr = VecDestroy(&localnearnullsp[k]);CHKERRQ(ierr); 882 } 883 ierr = PetscFree(localnearnullsp);CHKERRQ(ierr); 884 PetscFunctionReturn(0); 885 } 886 887 #undef __FUNCT__ 888 #define __FUNCT__ "PCBDDCAnalyzeInterface" 889 PetscErrorCode PCBDDCAnalyzeInterface(PC pc) 890 { 891 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 892 PC_IS *pcis = (PC_IS*)pc->data; 893 Mat_IS *matis = (Mat_IS*)pc->pmat->data; 894 PetscInt bs,ierr,i,vertex_size; 895 PetscViewer viewer=pcbddc->dbg_viewer; 896 897 PetscFunctionBegin; 898 /* Init local Graph struct */ 899 ierr = PCBDDCGraphInit(pcbddc->mat_graph,matis->mapping);CHKERRQ(ierr); 900 901 /* Check validity of the csr graph passed in by the user */ 902 if (pcbddc->mat_graph->nvtxs_csr != pcbddc->mat_graph->nvtxs) { 903 ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr); 904 } 905 /* Set default CSR adjacency of local dofs if not provided by the user with PCBDDCSetLocalAdjacencyGraph */ 906 if (!pcbddc->mat_graph->xadj || !pcbddc->mat_graph->adjncy) { 907 Mat mat_adj; 908 const PetscInt *xadj,*adjncy; 909 PetscBool flg_row=PETSC_TRUE; 910 911 ierr = MatConvert(matis->A,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr); 912 ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr); 913 if (!flg_row) { 914 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__); 915 } 916 ierr = PCBDDCSetLocalAdjacencyGraph(pc,i,xadj,adjncy,PETSC_COPY_VALUES);CHKERRQ(ierr); 917 ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr); 918 if (!flg_row) { 919 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 920 } 921 ierr = MatDestroy(&mat_adj);CHKERRQ(ierr); 922 } 923 924 /* Set default dofs' splitting if no information has been provided by the user with PCBDDCSetDofsSplitting */ 925 vertex_size = 1; 926 if (!pcbddc->n_ISForDofs) { 927 IS *custom_ISForDofs; 928 929 ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 930 ierr = PetscMalloc(bs*sizeof(IS),&custom_ISForDofs);CHKERRQ(ierr); 931 for (i=0;i<bs;i++) { 932 ierr = ISCreateStride(PETSC_COMM_SELF,pcis->n/bs,i,bs,&custom_ISForDofs[i]);CHKERRQ(ierr); 933 } 934 ierr = PCBDDCSetDofsSplitting(pc,bs,custom_ISForDofs);CHKERRQ(ierr); 935 /* remove my references to IS objects */ 936 for (i=0;i<bs;i++) { 937 ierr = ISDestroy(&custom_ISForDofs[i]);CHKERRQ(ierr); 938 } 939 ierr = PetscFree(custom_ISForDofs);CHKERRQ(ierr); 940 } else { /* mat block size as vertex size (used for elasticity) */ 941 ierr = MatGetBlockSize(matis->A,&vertex_size);CHKERRQ(ierr); 942 } 943 944 /* Setup of Graph */ 945 ierr = PCBDDCGraphSetUp(pcbddc->mat_graph,vertex_size,pcbddc->NeumannBoundaries,pcbddc->DirichletBoundaries,pcbddc->n_ISForDofs,pcbddc->ISForDofs,pcbddc->user_primal_vertices); 946 947 /* Graph's connected components analysis */ 948 ierr = PCBDDCGraphComputeConnectedComponents(pcbddc->mat_graph);CHKERRQ(ierr); 949 950 /* print some info to stdout */ 951 if (pcbddc->dbg_flag) { 952 ierr = PCBDDCGraphASCIIView(pcbddc->mat_graph,pcbddc->dbg_flag,viewer); 953 } 954 PetscFunctionReturn(0); 955 } 956 957 #undef __FUNCT__ 958 #define __FUNCT__ "PCBDDCGetPrimalVerticesLocalIdx" 959 PetscErrorCode PCBDDCGetPrimalVerticesLocalIdx(PC pc, PetscInt *n_vertices, PetscInt *vertices_idx[]) 960 { 961 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 962 PetscInt *vertices,*row_cmat_indices,n,i,size_of_constraint,local_primal_size; 963 PetscErrorCode ierr; 964 965 PetscFunctionBegin; 966 n = 0; 967 vertices = 0; 968 if (pcbddc->ConstraintMatrix) { 969 ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&i);CHKERRQ(ierr); 970 for (i=0;i<local_primal_size;i++) { 971 ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 972 if (size_of_constraint == 1) n++; 973 ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 974 } 975 ierr = PetscMalloc(n*sizeof(PetscInt),&vertices);CHKERRQ(ierr); 976 n = 0; 977 for (i=0;i<local_primal_size;i++) { 978 ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 979 if (size_of_constraint == 1) { 980 vertices[n++]=row_cmat_indices[0]; 981 } 982 ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 983 } 984 } 985 *n_vertices = n; 986 *vertices_idx = vertices; 987 PetscFunctionReturn(0); 988 } 989 990 #undef __FUNCT__ 991 #define __FUNCT__ "PCBDDCGetPrimalConstraintsLocalIdx" 992 PetscErrorCode PCBDDCGetPrimalConstraintsLocalIdx(PC pc, PetscInt *n_constraints, PetscInt *constraints_idx[]) 993 { 994 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 995 PetscInt *constraints_index,*row_cmat_indices,*row_cmat_global_indices; 996 PetscInt n,i,j,size_of_constraint,local_primal_size,local_size,max_size_of_constraint,min_index,min_loc; 997 PetscBool *touched; 998 PetscErrorCode ierr; 999 1000 PetscFunctionBegin; 1001 n = 0; 1002 constraints_index = 0; 1003 if (pcbddc->ConstraintMatrix) { 1004 ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&local_size);CHKERRQ(ierr); 1005 max_size_of_constraint = 0; 1006 for (i=0;i<local_primal_size;i++) { 1007 ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 1008 if (size_of_constraint > 1) { 1009 n++; 1010 } 1011 max_size_of_constraint = PetscMax(size_of_constraint,max_size_of_constraint); 1012 ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 1013 } 1014 ierr = PetscMalloc(n*sizeof(PetscInt),&constraints_index);CHKERRQ(ierr); 1015 ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&row_cmat_global_indices);CHKERRQ(ierr); 1016 ierr = PetscMalloc(local_size*sizeof(PetscBool),&touched);CHKERRQ(ierr); 1017 ierr = PetscMemzero(touched,local_size*sizeof(PetscBool));CHKERRQ(ierr); 1018 n = 0; 1019 for (i=0;i<local_primal_size;i++) { 1020 ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1021 if (size_of_constraint > 1) { 1022 ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,row_cmat_indices,row_cmat_global_indices);CHKERRQ(ierr); 1023 min_index = row_cmat_global_indices[0]; 1024 min_loc = 0; 1025 for (j=1;j<size_of_constraint;j++) { 1026 /* there can be more than one constraint on a single connected component */ 1027 if (min_index > row_cmat_global_indices[j] && !touched[row_cmat_indices[j]]) { 1028 min_index = row_cmat_global_indices[j]; 1029 min_loc = j; 1030 } 1031 } 1032 touched[row_cmat_indices[min_loc]] = PETSC_TRUE; 1033 constraints_index[n++] = row_cmat_indices[min_loc]; 1034 } 1035 ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1036 } 1037 } 1038 ierr = PetscFree(touched);CHKERRQ(ierr); 1039 ierr = PetscFree(row_cmat_global_indices);CHKERRQ(ierr); 1040 *n_constraints = n; 1041 *constraints_idx = constraints_index; 1042 PetscFunctionReturn(0); 1043 } 1044 1045 /* the next two functions has been adapted from pcis.c */ 1046 #undef __FUNCT__ 1047 #define __FUNCT__ "PCBDDCApplySchur" 1048 PetscErrorCode PCBDDCApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 1049 { 1050 PetscErrorCode ierr; 1051 PC_IS *pcis = (PC_IS*)(pc->data); 1052 1053 PetscFunctionBegin; 1054 if (!vec2_B) { vec2_B = v; } 1055 ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 1056 ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr); 1057 ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 1058 ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr); 1059 ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 1060 PetscFunctionReturn(0); 1061 } 1062 1063 #undef __FUNCT__ 1064 #define __FUNCT__ "PCBDDCApplySchurTranspose" 1065 PetscErrorCode PCBDDCApplySchurTranspose(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 1066 { 1067 PetscErrorCode ierr; 1068 PC_IS *pcis = (PC_IS*)(pc->data); 1069 1070 PetscFunctionBegin; 1071 if (!vec2_B) { vec2_B = v; } 1072 ierr = MatMultTranspose(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 1073 ierr = MatMultTranspose(pcis->A_BI,v,vec1_D);CHKERRQ(ierr); 1074 ierr = KSPSolveTranspose(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 1075 ierr = MatMultTranspose(pcis->A_IB,vec2_D,vec2_B);CHKERRQ(ierr); 1076 ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 1077 PetscFunctionReturn(0); 1078 } 1079 1080 #undef __FUNCT__ 1081 #define __FUNCT__ "PCBDDCSubsetNumbering" 1082 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[]) 1083 { 1084 Vec local_vec,global_vec; 1085 IS seqis,paris; 1086 VecScatter scatter_ctx; 1087 PetscScalar *array; 1088 PetscInt *temp_global_dofs; 1089 PetscScalar globalsum; 1090 PetscInt i,j,s; 1091 PetscInt nlocals,first_index,old_index,max_local; 1092 PetscMPIInt rank_prec_comm,size_prec_comm,max_global; 1093 PetscMPIInt *dof_sizes,*dof_displs; 1094 PetscBool first_found; 1095 PetscErrorCode ierr; 1096 1097 PetscFunctionBegin; 1098 /* mpi buffers */ 1099 MPI_Comm_size(comm,&size_prec_comm); 1100 MPI_Comm_rank(comm,&rank_prec_comm); 1101 j = ( !rank_prec_comm ? size_prec_comm : 0); 1102 ierr = PetscMalloc(j*sizeof(*dof_sizes),&dof_sizes);CHKERRQ(ierr); 1103 ierr = PetscMalloc(j*sizeof(*dof_displs),&dof_displs);CHKERRQ(ierr); 1104 /* get maximum size of subset */ 1105 ierr = PetscMalloc(n_local_dofs*sizeof(PetscInt),&temp_global_dofs);CHKERRQ(ierr); 1106 ierr = ISLocalToGlobalMappingApply(l2gmap,n_local_dofs,local_dofs,temp_global_dofs);CHKERRQ(ierr); 1107 max_local = 0; 1108 if (n_local_dofs) { 1109 max_local = temp_global_dofs[0]; 1110 for (i=1;i<n_local_dofs;i++) { 1111 if (max_local < temp_global_dofs[i] ) { 1112 max_local = temp_global_dofs[i]; 1113 } 1114 } 1115 } 1116 ierr = MPI_Allreduce(&max_local,&max_global,1,MPIU_INT,MPI_MAX,comm); 1117 max_global++; 1118 max_local = 0; 1119 if (n_local_dofs) { 1120 max_local = local_dofs[0]; 1121 for (i=1;i<n_local_dofs;i++) { 1122 if (max_local < local_dofs[i] ) { 1123 max_local = local_dofs[i]; 1124 } 1125 } 1126 } 1127 max_local++; 1128 /* allocate workspace */ 1129 ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr); 1130 ierr = VecSetSizes(local_vec,PETSC_DECIDE,max_local);CHKERRQ(ierr); 1131 ierr = VecSetType(local_vec,VECSEQ);CHKERRQ(ierr); 1132 ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr); 1133 ierr = VecSetSizes(global_vec,PETSC_DECIDE,max_global);CHKERRQ(ierr); 1134 ierr = VecSetType(global_vec,VECMPI);CHKERRQ(ierr); 1135 /* create scatter */ 1136 ierr = ISCreateGeneral(PETSC_COMM_SELF,n_local_dofs,local_dofs,PETSC_COPY_VALUES,&seqis);CHKERRQ(ierr); 1137 ierr = ISCreateGeneral(comm,n_local_dofs,temp_global_dofs,PETSC_COPY_VALUES,&paris);CHKERRQ(ierr); 1138 ierr = VecScatterCreate(local_vec,seqis,global_vec,paris,&scatter_ctx);CHKERRQ(ierr); 1139 ierr = ISDestroy(&seqis);CHKERRQ(ierr); 1140 ierr = ISDestroy(&paris);CHKERRQ(ierr); 1141 /* init array */ 1142 ierr = VecSet(global_vec,0.0);CHKERRQ(ierr); 1143 ierr = VecSet(local_vec,0.0);CHKERRQ(ierr); 1144 ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 1145 if (local_dofs_mult) { 1146 for (i=0;i<n_local_dofs;i++) { 1147 array[local_dofs[i]]=(PetscScalar)local_dofs_mult[i]; 1148 } 1149 } else { 1150 for (i=0;i<n_local_dofs;i++) { 1151 array[local_dofs[i]]=1.0; 1152 } 1153 } 1154 ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 1155 /* scatter into global vec and get total number of global dofs */ 1156 ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1157 ierr = VecScatterEnd (scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1158 ierr = VecSum(global_vec,&globalsum);CHKERRQ(ierr); 1159 *n_global_subset = (PetscInt)globalsum; 1160 /* Fill global_vec with cumulative function for global numbering */ 1161 ierr = VecGetArray(global_vec,&array);CHKERRQ(ierr); 1162 ierr = VecGetLocalSize(global_vec,&s);CHKERRQ(ierr); 1163 nlocals = 0; 1164 first_index = -1; 1165 first_found = PETSC_FALSE; 1166 for (i=0;i<s;i++) { 1167 if (!first_found && array[i] > 0.0) { 1168 first_found = PETSC_TRUE; 1169 first_index = i; 1170 } 1171 nlocals += (PetscInt)array[i]; 1172 } 1173 ierr = MPI_Gather(&nlocals,1,MPIU_INT,dof_sizes,1,MPIU_INT,0,comm);CHKERRQ(ierr); 1174 if (!rank_prec_comm) { 1175 dof_displs[0]=0; 1176 for (i=1;i<size_prec_comm;i++) { 1177 dof_displs[i] = dof_displs[i-1]+dof_sizes[i-1]; 1178 } 1179 } 1180 ierr = MPI_Scatter(dof_displs,1,MPIU_INT,&nlocals,1,MPIU_INT,0,comm);CHKERRQ(ierr); 1181 if (first_found) { 1182 array[first_index] += (PetscScalar)nlocals; 1183 old_index = first_index; 1184 for (i=first_index+1;i<s;i++) { 1185 if (array[i] > 0.0) { 1186 array[i] += array[old_index]; 1187 old_index = i; 1188 } 1189 } 1190 } 1191 ierr = VecRestoreArray(global_vec,&array);CHKERRQ(ierr); 1192 ierr = VecSet(local_vec,0.0);CHKERRQ(ierr); 1193 ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1194 ierr = VecScatterEnd (scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1195 /* get global ordering of local dofs */ 1196 ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 1197 if (local_dofs_mult) { 1198 for (i=0;i<n_local_dofs;i++) { 1199 temp_global_dofs[i] = (PetscInt)array[local_dofs[i]]-local_dofs_mult[i]; 1200 } 1201 } else { 1202 for (i=0;i<n_local_dofs;i++) { 1203 temp_global_dofs[i] = (PetscInt)array[local_dofs[i]]-1; 1204 } 1205 } 1206 ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 1207 /* free workspace */ 1208 ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr); 1209 ierr = VecDestroy(&local_vec);CHKERRQ(ierr); 1210 ierr = VecDestroy(&global_vec);CHKERRQ(ierr); 1211 ierr = PetscFree(dof_sizes);CHKERRQ(ierr); 1212 ierr = PetscFree(dof_displs);CHKERRQ(ierr); 1213 /* return pointer to global ordering of local dofs */ 1214 *global_numbering_subset = temp_global_dofs; 1215 PetscFunctionReturn(0); 1216 } 1217