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 /* uncomment for testing purposes */ 260 /* #define PETSC_MISSING_LAPACK_GESVD 1 */ 261 #undef __FUNCT__ 262 #define __FUNCT__ "PCBDDCConstraintsSetUp" 263 PetscErrorCode PCBDDCConstraintsSetUp(PC pc) 264 { 265 PetscErrorCode ierr; 266 PC_IS* pcis = (PC_IS*)(pc->data); 267 PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 268 Mat_IS* matis = (Mat_IS*)pc->pmat->data; 269 /* constraint and (optionally) change of basis matrix implemented as SeqAIJ */ 270 MatType impMatType=MATSEQAIJ; 271 /* one and zero */ 272 PetscScalar one=1.0,zero=0.0; 273 /* space to store constraints and their local indices */ 274 PetscScalar *temp_quadrature_constraint; 275 PetscInt *temp_indices,*temp_indices_to_constraint,*temp_indices_to_constraint_B; 276 /* iterators */ 277 PetscInt i,j,k,total_counts,temp_start_ptr; 278 /* stuff to store connected components stored in pcbddc->mat_graph */ 279 IS ISForVertices,*ISForFaces,*ISForEdges,*used_IS; 280 PetscInt n_ISForFaces,n_ISForEdges; 281 PetscBool get_faces,get_edges,get_vertices; 282 /* near null space stuff */ 283 MatNullSpace nearnullsp; 284 const Vec *nearnullvecs; 285 Vec *localnearnullsp; 286 PetscBool nnsp_has_cnst; 287 PetscInt nnsp_size; 288 PetscScalar *array; 289 /* BLAS integers */ 290 PetscBLASInt Bs,Bt,lwork,lierr,Bone=1; 291 /* LAPACK working arrays for SVD or POD */ 292 PetscBool skip_lapack; 293 PetscScalar *work; 294 PetscReal *singular_vals; 295 #if defined(PETSC_USE_COMPLEX) 296 PetscReal *rwork; 297 #endif 298 #if defined(PETSC_MISSING_LAPACK_GESVD) 299 PetscScalar *temp_basis,*correlation_mat; 300 #endif 301 /* change of basis */ 302 PetscInt *aux_primal_numbering,*aux_primal_minloc,*global_indices; 303 PetscBool boolforchange,*change_basis,*touched; 304 /* auxiliary stuff */ 305 PetscInt *nnz,*is_indices,*local_to_B; 306 /* some quantities */ 307 PetscInt n_vertices,total_primal_vertices; 308 PetscInt size_of_constraint,max_size_of_constraint,max_constraints,temp_constraints; 309 310 311 PetscFunctionBegin; 312 /* Get index sets for faces, edges and vertices from graph */ 313 get_faces = PETSC_TRUE; 314 get_edges = PETSC_TRUE; 315 get_vertices = PETSC_TRUE; 316 if (pcbddc->vertices_flag) { 317 get_faces = PETSC_FALSE; 318 get_edges = PETSC_FALSE; 319 } 320 if (pcbddc->constraints_flag) { 321 get_vertices = PETSC_FALSE; 322 } 323 if (pcbddc->faces_flag) { 324 get_edges = PETSC_FALSE; 325 } 326 if (pcbddc->edges_flag) { 327 get_faces = PETSC_FALSE; 328 } 329 /* default */ 330 if (!get_faces && !get_edges && !get_vertices) { 331 get_vertices = PETSC_TRUE; 332 } 333 ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,get_faces,get_edges,get_vertices,&n_ISForFaces,&ISForFaces,&n_ISForEdges,&ISForEdges,&ISForVertices); 334 /* print some info */ 335 if (pcbddc->dbg_flag) { 336 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr); 337 i = 0; 338 if (ISForVertices) { 339 ierr = ISGetSize(ISForVertices,&i);CHKERRQ(ierr); 340 } 341 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate vertices\n",PetscGlobalRank,i);CHKERRQ(ierr); 342 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate edges\n",PetscGlobalRank,n_ISForEdges);CHKERRQ(ierr); 343 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate faces\n",PetscGlobalRank,n_ISForFaces);CHKERRQ(ierr); 344 ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 345 } 346 /* check if near null space is attached to global mat */ 347 ierr = MatGetNearNullSpace(pc->pmat,&nearnullsp);CHKERRQ(ierr); 348 if (nearnullsp) { 349 ierr = MatNullSpaceGetVecs(nearnullsp,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr); 350 } else { /* if near null space is not provided BDDC uses constants by default */ 351 nnsp_size = 0; 352 nnsp_has_cnst = PETSC_TRUE; 353 } 354 /* get max number of constraints on a single cc */ 355 max_constraints = nnsp_size; 356 if (nnsp_has_cnst) max_constraints++; 357 358 /* 359 Evaluate maximum storage size needed by the procedure 360 - temp_indices will contain start index of each constraint stored as follows 361 - temp_indices_to_constraint [temp_indices[i],...,temp[indices[i+1]-1] will contain the indices (in local numbering) on which the constraint acts 362 - 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 363 - temp_quadrature_constraint [temp_indices[i],...,temp[indices[i+1]-1] will contain the scalars representing the constraint itself 364 */ 365 total_counts = n_ISForFaces+n_ISForEdges; 366 total_counts *= max_constraints; 367 n_vertices = 0; 368 if (ISForVertices) { 369 ierr = ISGetSize(ISForVertices,&n_vertices);CHKERRQ(ierr); 370 } 371 total_counts += n_vertices; 372 ierr = PetscMalloc((total_counts+1)*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr); 373 ierr = PetscMalloc((total_counts+1)*sizeof(PetscBool),&change_basis);CHKERRQ(ierr); 374 total_counts = 0; 375 max_size_of_constraint = 0; 376 for (i=0;i<n_ISForEdges+n_ISForFaces;i++) { 377 if (i<n_ISForEdges) { 378 used_IS = &ISForEdges[i]; 379 } else { 380 used_IS = &ISForFaces[i-n_ISForEdges]; 381 } 382 ierr = ISGetSize(*used_IS,&j);CHKERRQ(ierr); 383 total_counts += j; 384 max_size_of_constraint = PetscMax(j,max_size_of_constraint); 385 } 386 total_counts *= max_constraints; 387 total_counts += n_vertices; 388 ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&temp_quadrature_constraint);CHKERRQ(ierr); 389 ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint);CHKERRQ(ierr); 390 ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint_B);CHKERRQ(ierr); 391 /* local to boundary numbering */ 392 ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&local_to_B);CHKERRQ(ierr); 393 ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 394 for (i=0;i<pcis->n;i++) local_to_B[i]=-1; 395 for (i=0;i<pcis->n_B;i++) local_to_B[is_indices[i]]=i; 396 ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 397 /* get local part of global near null space vectors */ 398 ierr = PetscMalloc(nnsp_size*sizeof(Vec),&localnearnullsp);CHKERRQ(ierr); 399 for (k=0;k<nnsp_size;k++) { 400 ierr = VecDuplicate(pcis->vec1_N,&localnearnullsp[k]);CHKERRQ(ierr); 401 ierr = VecScatterBegin(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 402 ierr = VecScatterEnd(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 403 } 404 405 /* whether or not to skip lapack calls */ 406 skip_lapack = PETSC_TRUE; 407 if (n_ISForFaces+n_ISForEdges) skip_lapack = PETSC_FALSE; 408 409 /* First we issue queries to allocate optimal workspace for LAPACKgesvd (or LAPACKsyev if SVD is missing) */ 410 if (!pcbddc->use_nnsp_true && !skip_lapack) { 411 PetscScalar temp_work; 412 #if defined(PETSC_MISSING_LAPACK_GESVD) 413 /* Proper Orthogonal Decomposition (POD) using the snapshot method */ 414 ierr = PetscMalloc(max_constraints*max_constraints*sizeof(PetscScalar),&correlation_mat);CHKERRQ(ierr); 415 ierr = PetscMalloc(max_constraints*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr); 416 ierr = PetscMalloc(max_size_of_constraint*max_constraints*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr); 417 #if defined(PETSC_USE_COMPLEX) 418 ierr = PetscMalloc(3*max_constraints*sizeof(PetscReal),&rwork);CHKERRQ(ierr); 419 #endif 420 /* now we evaluate the optimal workspace using query with lwork=-1 */ 421 ierr = PetscBLASIntCast(max_constraints,&Bt);CHKERRQ(ierr); 422 lwork = -1; 423 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 424 #if !defined(PETSC_USE_COMPLEX) 425 PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,&lierr)); 426 #else 427 PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,rwork,&lierr)); 428 #endif 429 ierr = PetscFPTrapPop();CHKERRQ(ierr); 430 if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYEV Lapack routine %d",(int)lierr); 431 #else /* on missing GESVD */ 432 /* SVD */ 433 PetscInt max_n,min_n; 434 max_n = max_size_of_constraint; 435 min_n = max_constraints; 436 if (max_size_of_constraint < max_constraints) { 437 min_n = max_size_of_constraint; 438 max_n = max_constraints; 439 } 440 ierr = PetscMalloc(min_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr); 441 #if defined(PETSC_USE_COMPLEX) 442 ierr = PetscMalloc(5*min_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr); 443 #endif 444 /* now we evaluate the optimal workspace using query with lwork=-1 */ 445 lwork = -1; 446 ierr = PetscBLASIntCast(max_n,&Bs);CHKERRQ(ierr); 447 ierr = PetscBLASIntCast(min_n,&Bt);CHKERRQ(ierr); 448 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 449 #if !defined(PETSC_USE_COMPLEX) 450 PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Bs,singular_vals,work,&Bt,work,&Bt,&temp_work,&lwork,&lierr)); 451 #else 452 PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Bs,singular_vals,work,&Bt,work,&Bt,&temp_work,&lwork,rwork,&lierr)); 453 #endif 454 ierr = PetscFPTrapPop();CHKERRQ(ierr); 455 if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GESVD Lapack routine %d",(int)lierr); 456 #endif /* on missing GESVD */ 457 /* Allocate optimal workspace */ 458 ierr = PetscBLASIntCast((PetscInt)PetscRealPart(temp_work),&lwork);CHKERRQ(ierr); 459 ierr = PetscMalloc((PetscInt)lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr); 460 } 461 /* Now we can loop on constraining sets */ 462 total_counts = 0; 463 temp_indices[0] = 0; 464 /* vertices */ 465 if (ISForVertices) { 466 ierr = ISGetIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); 467 if (nnsp_has_cnst) { /* consider all vertices */ 468 for (i=0;i<n_vertices;i++) { 469 temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i]; 470 temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]]; 471 temp_quadrature_constraint[temp_indices[total_counts]]=1.0; 472 temp_indices[total_counts+1]=temp_indices[total_counts]+1; 473 change_basis[total_counts]=PETSC_FALSE; 474 total_counts++; 475 } 476 } else { /* consider vertices for which exist at least a localnearnullsp which is not null there */ 477 PetscBool used_vertex; 478 for (i=0;i<n_vertices;i++) { 479 used_vertex = PETSC_FALSE; 480 k = 0; 481 while (!used_vertex && k<nnsp_size) { 482 ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr); 483 if (PetscAbsScalar(array[is_indices[i]])>0.0) { 484 temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i]; 485 temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]]; 486 temp_quadrature_constraint[temp_indices[total_counts]]=1.0; 487 temp_indices[total_counts+1]=temp_indices[total_counts]+1; 488 change_basis[total_counts]=PETSC_FALSE; 489 total_counts++; 490 used_vertex = PETSC_TRUE; 491 } 492 ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr); 493 k++; 494 } 495 } 496 } 497 ierr = ISRestoreIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); 498 n_vertices = total_counts; 499 } 500 501 /* edges and faces */ 502 for (i=0;i<n_ISForEdges+n_ISForFaces;i++) { 503 if (i<n_ISForEdges) { 504 used_IS = &ISForEdges[i]; 505 boolforchange = pcbddc->use_change_of_basis; /* change or not the basis on the edge */ 506 } else { 507 used_IS = &ISForFaces[i-n_ISForEdges]; 508 boolforchange = (PetscBool)(pcbddc->use_change_of_basis && pcbddc->use_change_on_faces); /* change or not the basis on the face */ 509 } 510 temp_constraints = 0; /* zero the number of constraints I have on this conn comp */ 511 temp_start_ptr = total_counts; /* need to know the starting index of constraints stored */ 512 ierr = ISGetSize(*used_IS,&size_of_constraint);CHKERRQ(ierr); 513 ierr = ISGetIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 514 /* change of basis should not be performed on local periodic nodes */ 515 if (pcbddc->mat_graph->mirrors && pcbddc->mat_graph->mirrors[is_indices[0]]) boolforchange = PETSC_FALSE; 516 if (nnsp_has_cnst) { 517 PetscScalar quad_value; 518 temp_constraints++; 519 quad_value = (PetscScalar)(1.0/PetscSqrtReal((PetscReal)size_of_constraint)); 520 for (j=0;j<size_of_constraint;j++) { 521 temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j]; 522 temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]]; 523 temp_quadrature_constraint[temp_indices[total_counts]+j]=quad_value; 524 } 525 temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint; /* store new starting point */ 526 change_basis[total_counts]=boolforchange; 527 total_counts++; 528 } 529 for (k=0;k<nnsp_size;k++) { 530 PetscReal real_value; 531 ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr); 532 for (j=0;j<size_of_constraint;j++) { 533 temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j]; 534 temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]]; 535 temp_quadrature_constraint[temp_indices[total_counts]+j]=array[is_indices[j]]; 536 } 537 ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr); 538 /* check if array is null on the connected component */ 539 ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr); 540 PetscStackCallBLAS("BLASasum",real_value = BLASasum_(&Bs,&temp_quadrature_constraint[temp_indices[total_counts]],&Bone)); 541 if (real_value > 0.0) { /* keep indices and values */ 542 temp_constraints++; 543 temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint; /* store new starting point */ 544 change_basis[total_counts]=boolforchange; 545 total_counts++; 546 } 547 } 548 ierr = ISRestoreIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 549 /* perform SVD on the constraints if use_nnsp_true has not be requested by the user */ 550 if (!pcbddc->use_nnsp_true) { 551 PetscReal tol = 1.0e-8; /* tolerance for retaining eigenmodes */ 552 553 ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr); 554 ierr = PetscBLASIntCast(temp_constraints,&Bt);CHKERRQ(ierr); 555 #if defined(PETSC_MISSING_LAPACK_GESVD) 556 /* SVD: Y = U*S*V^H -> U (eigenvectors of Y*Y^H) = Y*V*(S)^\dag 557 POD: Y^H*Y = V*D*V^H, D = S^H*S -> U = Y*V*D^(-1/2) 558 -> When PETSC_USE_COMPLEX and PETSC_MISSING_LAPACK_GESVD are defined 559 the constraints basis will differ (by a complex factor with absolute value equal to 1) 560 from that computed using LAPACKgesvd 561 -> This is due to a different computation of eigenvectors in LAPACKheev 562 -> The quality of the POD-computed basis will be the same */ 563 ierr = PetscMemzero(correlation_mat,temp_constraints*temp_constraints*sizeof(PetscScalar));CHKERRQ(ierr); 564 /* Store upper triangular part of correlation matrix */ 565 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 566 for (j=0;j<temp_constraints;j++) { 567 for (k=0;k<j+1;k++) { 568 PetscStackCallBLAS("BLASdot",correlation_mat[j*temp_constraints+k]=BLASdot_(&Bs,&temp_quadrature_constraint[temp_indices[temp_start_ptr+k]],&Bone,&temp_quadrature_constraint[temp_indices[temp_start_ptr+j]],&Bone)); 569 } 570 } 571 #if !defined(PETSC_USE_COMPLEX) 572 PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,work,&lwork,&lierr)); 573 #else 574 PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,work,&lwork,rwork,&lierr)); 575 #endif 576 ierr = PetscFPTrapPop();CHKERRQ(ierr); 577 if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYEV Lapack routine %d",(int)lierr); 578 /* retain eigenvalues greater than tol: note that LAPACKsyev gives eigs in ascending order */ 579 j=0; 580 while (j < temp_constraints && singular_vals[j] < tol) j++; 581 total_counts=total_counts-j; 582 if (j<temp_constraints) { 583 PetscInt ii; 584 for (k=j;k<temp_constraints;k++) singular_vals[k]=1.0/PetscSqrtReal(singular_vals[k]); 585 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 586 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)); 587 ierr = PetscFPTrapPop();CHKERRQ(ierr); 588 /* scale and copy POD basis into used quadrature memory */ 589 for (k=0;k<temp_constraints-j;k++) { 590 for (ii=0;ii<size_of_constraint;ii++) { 591 temp_quadrature_constraint[temp_indices[temp_start_ptr+k]+ii]=singular_vals[temp_constraints-1-k]*temp_basis[(temp_constraints-1-k)*size_of_constraint+ii]; 592 } 593 } 594 } 595 #else /* on missing GESVD */ 596 PetscInt min_n = temp_constraints; 597 if (min_n > size_of_constraint) min_n = size_of_constraint; 598 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 599 #if !defined(PETSC_USE_COMPLEX) 600 PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,singular_vals,work,&Bt,work,&Bt,work,&lwork,&lierr)); 601 #else 602 PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,singular_vals,work,&Bt,work,&Bt,work,&lwork,rwork,&lierr)); 603 #endif 604 if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GESVD Lapack routine %d",(int)lierr); 605 ierr = PetscFPTrapPop();CHKERRQ(ierr); 606 /* retain eigenvalues greater than tol: note that LAPACKgesvd gives eigs in descending order */ 607 j = 0; 608 while (j < min_n && singular_vals[min_n-j-1] < tol) j++; 609 total_counts = total_counts-temp_constraints+min_n-j; 610 #endif /* on missing GESVD */ 611 } 612 } 613 /* free index sets of faces, edges and vertices */ 614 for (i=0;i<n_ISForFaces;i++) { 615 ierr = ISDestroy(&ISForFaces[i]);CHKERRQ(ierr); 616 } 617 ierr = PetscFree(ISForFaces);CHKERRQ(ierr); 618 for (i=0;i<n_ISForEdges;i++) { 619 ierr = ISDestroy(&ISForEdges[i]);CHKERRQ(ierr); 620 } 621 ierr = PetscFree(ISForEdges);CHKERRQ(ierr); 622 ierr = ISDestroy(&ISForVertices);CHKERRQ(ierr); 623 624 /* free workspace */ 625 if (!pcbddc->use_nnsp_true && !skip_lapack) { 626 ierr = PetscFree(work);CHKERRQ(ierr); 627 #if defined(PETSC_USE_COMPLEX) 628 ierr = PetscFree(rwork);CHKERRQ(ierr); 629 #endif 630 ierr = PetscFree(singular_vals);CHKERRQ(ierr); 631 #if defined(PETSC_MISSING_LAPACK_GESVD) 632 ierr = PetscFree(correlation_mat);CHKERRQ(ierr); 633 ierr = PetscFree(temp_basis);CHKERRQ(ierr); 634 #endif 635 } 636 for (k=0;k<nnsp_size;k++) { 637 ierr = VecDestroy(&localnearnullsp[k]);CHKERRQ(ierr); 638 } 639 ierr = PetscFree(localnearnullsp);CHKERRQ(ierr); 640 641 /* set quantities in pcbddc data structure */ 642 /* n_vertices defines the number of subdomain corners in the primal space */ 643 /* n_constraints defines the number of averages (they can be point primal dofs if change of basis is requested) */ 644 pcbddc->local_primal_size = total_counts; 645 pcbddc->n_vertices = n_vertices; 646 pcbddc->n_constraints = pcbddc->local_primal_size-pcbddc->n_vertices; 647 648 /* Create constraint matrix */ 649 /* The constraint matrix is used to compute the l2g map of primal dofs */ 650 /* so we need to set it up properly either with or without change of basis */ 651 ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ConstraintMatrix);CHKERRQ(ierr); 652 ierr = MatSetType(pcbddc->ConstraintMatrix,impMatType);CHKERRQ(ierr); 653 ierr = MatSetSizes(pcbddc->ConstraintMatrix,pcbddc->local_primal_size,pcis->n,pcbddc->local_primal_size,pcis->n);CHKERRQ(ierr); 654 /* array to compute a local numbering of constraints : vertices first then constraints */ 655 ierr = PetscMalloc(pcbddc->local_primal_size*sizeof(PetscInt),&aux_primal_numbering);CHKERRQ(ierr); 656 /* array to select the proper local node (of minimum index with respect to global ordering) when changing the basis */ 657 /* note: it should not be needed since IS for faces and edges are already sorted by global ordering when analyzing the graph but... just in case */ 658 ierr = PetscMalloc(pcbddc->local_primal_size*sizeof(PetscInt),&aux_primal_minloc);CHKERRQ(ierr); 659 /* auxiliary stuff for basis change */ 660 ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&global_indices);CHKERRQ(ierr); 661 ierr = PetscMalloc(pcis->n_B*sizeof(PetscBool),&touched);CHKERRQ(ierr); 662 ierr = PetscMemzero(touched,pcis->n_B*sizeof(PetscBool));CHKERRQ(ierr); 663 664 /* find primal_dofs: subdomain corners plus dofs selected as primal after change of basis */ 665 total_primal_vertices=0; 666 for (i=0;i<pcbddc->local_primal_size;i++) { 667 size_of_constraint=temp_indices[i+1]-temp_indices[i]; 668 if (size_of_constraint == 1) { 669 touched[temp_indices_to_constraint_B[temp_indices[i]]]=PETSC_TRUE; 670 aux_primal_numbering[total_primal_vertices]=temp_indices_to_constraint[temp_indices[i]]; 671 aux_primal_minloc[total_primal_vertices]=0; 672 total_primal_vertices++; 673 } else if (change_basis[i]) { /* Same procedure used in PCBDDCGetPrimalConstraintsLocalIdx */ 674 PetscInt min_loc,min_index; 675 ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,&temp_indices_to_constraint[temp_indices[i]],global_indices);CHKERRQ(ierr); 676 /* find first untouched local node */ 677 k = 0; 678 while (touched[temp_indices_to_constraint_B[temp_indices[i]+k]]) k++; 679 min_index = global_indices[k]; 680 min_loc = k; 681 /* search the minimum among global nodes already untouched on the cc */ 682 for (k=1;k<size_of_constraint;k++) { 683 /* there can be more than one constraint on a single connected component */ 684 if (min_index > global_indices[k] && !touched[temp_indices_to_constraint_B[temp_indices[i]+k]]) { 685 min_index = global_indices[k]; 686 min_loc = k; 687 } 688 } 689 touched[temp_indices_to_constraint_B[temp_indices[i]+min_loc]] = PETSC_TRUE; 690 aux_primal_numbering[total_primal_vertices]=temp_indices_to_constraint[temp_indices[i]+min_loc]; 691 aux_primal_minloc[total_primal_vertices]=min_loc; 692 total_primal_vertices++; 693 } 694 } 695 /* free workspace */ 696 ierr = PetscFree(global_indices);CHKERRQ(ierr); 697 ierr = PetscFree(touched);CHKERRQ(ierr); 698 /* permute indices in order to have a sorted set of vertices */ 699 ierr = PetscSortInt(total_primal_vertices,aux_primal_numbering); 700 701 /* nonzero structure of constraint matrix */ 702 ierr = PetscMalloc(pcbddc->local_primal_size*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 703 for (i=0;i<total_primal_vertices;i++) nnz[i]=1; 704 j=total_primal_vertices; 705 for (i=pcbddc->n_vertices;i<pcbddc->local_primal_size;i++) { 706 if (!change_basis[i]) { 707 nnz[j]=temp_indices[i+1]-temp_indices[i]; 708 j++; 709 } 710 } 711 ierr = MatSeqAIJSetPreallocation(pcbddc->ConstraintMatrix,0,nnz);CHKERRQ(ierr); 712 ierr = PetscFree(nnz);CHKERRQ(ierr); 713 /* set values in constraint matrix */ 714 for (i=0;i<total_primal_vertices;i++) { 715 ierr = MatSetValue(pcbddc->ConstraintMatrix,i,aux_primal_numbering[i],1.0,INSERT_VALUES);CHKERRQ(ierr); 716 } 717 total_counts = total_primal_vertices; 718 for (i=pcbddc->n_vertices;i<pcbddc->local_primal_size;i++) { 719 if (!change_basis[i]) { 720 size_of_constraint=temp_indices[i+1]-temp_indices[i]; 721 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); 722 total_counts++; 723 } 724 } 725 /* assembling */ 726 ierr = MatAssemblyBegin(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 727 ierr = MatAssemblyEnd(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 728 /* 729 ierr = MatView(pcbddc->ConstraintMatrix,(PetscViewer)0);CHKERRQ(ierr); 730 */ 731 /* Create matrix for change of basis. We don't need it in case pcbddc->use_change_of_basis is FALSE */ 732 if (pcbddc->use_change_of_basis) { 733 PetscBool qr_needed = PETSC_FALSE; 734 /* change of basis acts on local interfaces -> dimension is n_B x n_B */ 735 ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr); 736 ierr = MatSetType(pcbddc->ChangeOfBasisMatrix,impMatType);CHKERRQ(ierr); 737 ierr = MatSetSizes(pcbddc->ChangeOfBasisMatrix,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);CHKERRQ(ierr); 738 /* work arrays */ 739 ierr = PetscMalloc(pcis->n_B*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 740 for (i=0;i<pcis->n_B;i++) nnz[i]=1; 741 /* nonzeros per row */ 742 for (i=pcbddc->n_vertices;i<pcbddc->local_primal_size;i++) { 743 if (change_basis[i]) { 744 qr_needed = PETSC_TRUE; 745 size_of_constraint = temp_indices[i+1]-temp_indices[i]; 746 for (j=0;j<size_of_constraint;j++) nnz[temp_indices_to_constraint_B[temp_indices[i]+j]] = size_of_constraint; 747 } 748 } 749 ierr = MatSeqAIJSetPreallocation(pcbddc->ChangeOfBasisMatrix,0,nnz);CHKERRQ(ierr); 750 ierr = PetscFree(nnz);CHKERRQ(ierr); 751 /* Set initial identity in the matrix */ 752 for (i=0;i<pcis->n_B;i++) { 753 ierr = MatSetValue(pcbddc->ChangeOfBasisMatrix,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr); 754 } 755 756 /* Now we loop on the constraints which need a change of basis */ 757 /* Change of basis matrix is evaluated as the FIRST APPROACH in */ 758 /* Klawonn and Widlund, Dual-primal FETI-DP methods for linear elasticity, (see Sect 6.2.1) */ 759 /* Change of basis matrix T computed via QR decomposition of constraints */ 760 if (qr_needed) { 761 /* dual and primal dofs on a single cc */ 762 PetscInt dual_dofs,primal_dofs; 763 /* iterator on aux_primal_minloc (ordered as read from nearnullspace: vertices, edges and then constraints) */ 764 PetscInt primal_counter; 765 /* working stuff for GEQRF */ 766 PetscScalar *qr_basis,*qr_tau,*qr_work,lqr_work_t; 767 PetscBLASInt lqr_work; 768 /* working stuff for UNGQR */ 769 PetscScalar *gqr_work,lgqr_work_t; 770 PetscBLASInt lgqr_work; 771 /* working stuff for TRTRS */ 772 PetscScalar *trs_rhs; 773 /* pointers for values insertion into change of basis matrix */ 774 PetscInt *start_rows,*start_cols; 775 PetscScalar *start_vals; 776 /* working stuff for values insertion */ 777 PetscBool *is_primal; 778 779 /* space to store Q */ 780 ierr = PetscMalloc((max_size_of_constraint)*(max_size_of_constraint)*sizeof(PetscScalar),&qr_basis);CHKERRQ(ierr); 781 /* first we issue queries for optimal work */ 782 ierr = PetscBLASIntCast(max_size_of_constraint,&Bs);CHKERRQ(ierr); 783 ierr = PetscBLASIntCast(max_constraints,&Bt);CHKERRQ(ierr); 784 lqr_work = -1; 785 PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Bs,&Bt,qr_basis,&Bs,qr_tau,&lqr_work_t,&lqr_work,&lierr)); 786 if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GEQRF Lapack routine %d",(int)lierr); 787 ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lqr_work_t),&lqr_work);CHKERRQ(ierr); 788 ierr = PetscMalloc((PetscInt)PetscRealPart(lqr_work_t)*sizeof(*qr_work),&qr_work);CHKERRQ(ierr); 789 lgqr_work = -1; 790 if (Bt>Bs) Bt=Bs; /* adjust Bt just for computing optimal work */ 791 PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Bs,&Bs,&Bt,qr_basis,&Bs,qr_tau,&lgqr_work_t,&lgqr_work,&lierr)); 792 if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to UNGQR Lapack routine %d",(int)lierr); 793 ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lgqr_work_t),&lgqr_work);CHKERRQ(ierr); 794 ierr = PetscMalloc((PetscInt)PetscRealPart(lgqr_work_t)*sizeof(*gqr_work),&gqr_work);CHKERRQ(ierr); 795 /* array to store scaling factors for reflectors */ 796 ierr = PetscMalloc(max_constraints*sizeof(*qr_tau),&qr_tau);CHKERRQ(ierr); 797 /* array to store rhs and solution of triangular solver */ 798 ierr = PetscMalloc(max_constraints*max_constraints*sizeof(*trs_rhs),&trs_rhs);CHKERRQ(ierr); 799 /* array to store whether a node is primal or not */ 800 ierr = PetscMalloc(pcis->n_B*sizeof(*is_primal),&is_primal);CHKERRQ(ierr); 801 ierr = PetscMemzero(is_primal,pcis->n_B*sizeof(*is_primal));CHKERRQ(ierr); 802 for (i=0;i<total_primal_vertices;i++) is_primal[local_to_B[aux_primal_numbering[i]]] = PETSC_TRUE; 803 804 /* allocating workspace for check */ 805 if (pcbddc->dbg_flag) { 806 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr); 807 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Checking change of basis computation for subdomain %04d\n",PetscGlobalRank);CHKERRQ(ierr); 808 ierr = PetscMalloc(max_size_of_constraint*(max_constraints+max_size_of_constraint)*sizeof(*work),&work);CHKERRQ(ierr); 809 } 810 811 /* loop on constraints and see whether or not they need a change of basis */ 812 /* -> using implicit ordering contained in temp_indices data */ 813 total_counts = pcbddc->n_vertices; 814 primal_counter = total_counts; 815 while (total_counts<pcbddc->local_primal_size) { 816 primal_dofs = 1; 817 if (change_basis[total_counts]) { 818 /* get all constraints with same support: if more then one constraint is present on the cc then surely indices are stored contiguosly */ 819 while (total_counts+primal_dofs < pcbddc->local_primal_size && temp_indices_to_constraint_B[temp_indices[total_counts]] == temp_indices_to_constraint_B[temp_indices[total_counts+primal_dofs]]) { 820 primal_dofs++; 821 } 822 /* get constraint info */ 823 size_of_constraint = temp_indices[total_counts+1]-temp_indices[total_counts]; 824 dual_dofs = size_of_constraint-primal_dofs; 825 /* get BLAS dims */ 826 ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr); 827 ierr = PetscBLASIntCast(primal_dofs,&Bt);CHKERRQ(ierr); 828 829 /* copy quadrature constraints for change of basis check */ 830 if (pcbddc->dbg_flag) { 831 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Constraint %d to %d need a change of basis (size %d)\n",total_counts,total_counts+primal_dofs,size_of_constraint);CHKERRQ(ierr); 832 ierr = PetscMemcpy(work,&temp_quadrature_constraint[temp_indices[total_counts]],size_of_constraint*primal_dofs*sizeof(PetscScalar));CHKERRQ(ierr); 833 } 834 835 /* copy temporary constraints into larger work vector (in order to store all columns of Q) */ 836 ierr = PetscMemcpy(qr_basis,&temp_quadrature_constraint[temp_indices[total_counts]],size_of_constraint*primal_dofs*sizeof(PetscScalar));CHKERRQ(ierr); 837 838 /* compute QR decomposition of constraints */ 839 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 840 PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Bs,&Bt,qr_basis,&Bs,qr_tau,qr_work,&lqr_work,&lierr)); 841 if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GEQRF Lapack routine %d",(int)lierr); 842 ierr = PetscFPTrapPop();CHKERRQ(ierr); 843 844 /* explictly compute R^-T */ 845 ierr = PetscMemzero(trs_rhs,primal_dofs*primal_dofs*sizeof(*trs_rhs));CHKERRQ(ierr); 846 for (j=0;j<primal_dofs;j++) trs_rhs[j*(primal_dofs+1)] = 1.0; 847 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 848 PetscStackCallBLAS("LAPACKtrtrs",LAPACKtrtrs_("U","T","N",&Bt,&Bt,qr_basis,&Bs,trs_rhs,&Bt,&lierr)); 849 if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in TRTRS Lapack routine %d",(int)lierr); 850 ierr = PetscFPTrapPop();CHKERRQ(ierr); 851 852 /* explcitly compute all columns of Q (Q = [Q1 | Q2] ) overwriting QR factorization in qr_basis */ 853 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 854 PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Bs,&Bs,&Bt,qr_basis,&Bs,qr_tau,gqr_work,&lgqr_work,&lierr)); 855 if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in UNGQR Lapack routine %d",(int)lierr); 856 ierr = PetscFPTrapPop();CHKERRQ(ierr); 857 858 /* first primal_dofs columns of Q need to be re-scaled in order to be unitary w.r.t constraints 859 i.e. C_{pxn}*Q_{nxn} should be equal to [I_pxp | 0_pxd] (see check below) 860 where n=size_of_constraint, p=primal_dofs, d=dual_dofs (n=p+d), I and 0 identity and null matrix resp. */ 861 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 862 PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&Bs,&Bt,&Bt,&one,qr_basis,&Bs,trs_rhs,&Bt,&zero,&temp_quadrature_constraint[temp_indices[total_counts]],&Bs)); 863 ierr = PetscFPTrapPop();CHKERRQ(ierr); 864 ierr = PetscMemcpy(qr_basis,&temp_quadrature_constraint[temp_indices[total_counts]],size_of_constraint*primal_dofs*sizeof(PetscScalar));CHKERRQ(ierr); 865 866 /* insert values in change of basis matrix respecting global ordering of new primal dofs */ 867 start_rows = &temp_indices_to_constraint_B[temp_indices[total_counts]]; 868 /* insert cols for primal dofs */ 869 for (j=0;j<primal_dofs;j++) { 870 start_vals = &qr_basis[j*size_of_constraint]; 871 start_cols = &temp_indices_to_constraint_B[temp_indices[total_counts]+aux_primal_minloc[primal_counter+j]]; 872 ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,start_rows,1,start_cols,start_vals,INSERT_VALUES);CHKERRQ(ierr); 873 } 874 /* insert cols for dual dofs */ 875 for (j=0,k=0;j<dual_dofs;k++) { 876 if (!is_primal[temp_indices_to_constraint_B[temp_indices[total_counts]+k]]) { 877 start_vals = &qr_basis[(primal_dofs+j)*size_of_constraint]; 878 start_cols = &temp_indices_to_constraint_B[temp_indices[total_counts]+k]; 879 ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,start_rows,1,start_cols,start_vals,INSERT_VALUES);CHKERRQ(ierr); 880 j++; 881 } 882 } 883 884 /* check change of basis */ 885 if (pcbddc->dbg_flag) { 886 PetscInt ii,jj; 887 PetscBool valid_qr=PETSC_TRUE; 888 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 889 PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&Bt,&Bs,&Bs,&one,work,&Bs,qr_basis,&Bs,&zero,&work[size_of_constraint*primal_dofs],&Bt)); 890 ierr = PetscFPTrapPop();CHKERRQ(ierr); 891 for (jj=0;jj<size_of_constraint;jj++) { 892 for (ii=0;ii<primal_dofs;ii++) { 893 if (ii != jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]) > 1.e-12) valid_qr = PETSC_FALSE; 894 if (ii == jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]-1.0) > 1.e-12) valid_qr = PETSC_FALSE; 895 } 896 } 897 if (!valid_qr) { 898 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"\t-> wrong change of basis!\n",PetscGlobalRank);CHKERRQ(ierr); 899 for (jj=0;jj<size_of_constraint;jj++) { 900 for (ii=0;ii<primal_dofs;ii++) { 901 if (ii != jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]) > 1.e-12) { 902 PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"\tQr basis function %d is not orthogonal to constraint %d (%1.14e)!\n",jj,ii,PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii])); 903 } 904 if (ii == jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]-1.0) > 1.e-12) { 905 PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"\tQr basis function %d is not unitary w.r.t constraint %d (%1.14e)!\n",jj,ii,PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii])); 906 } 907 } 908 } 909 } else { 910 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"\t-> right change of basis!\n",PetscGlobalRank);CHKERRQ(ierr); 911 } 912 } 913 /* increment primal counter */ 914 primal_counter += primal_dofs; 915 } else { 916 if (pcbddc->dbg_flag) { 917 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Constraint %d does not need a change of basis (size %d)\n",total_counts,temp_indices[total_counts+1]-temp_indices[total_counts]);CHKERRQ(ierr); 918 } 919 } 920 /* increment constraint counter total_counts */ 921 total_counts += primal_dofs; 922 } 923 if (pcbddc->dbg_flag) { 924 ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 925 ierr = PetscFree(work);CHKERRQ(ierr); 926 } 927 /* free workspace */ 928 ierr = PetscFree(trs_rhs);CHKERRQ(ierr); 929 ierr = PetscFree(qr_tau);CHKERRQ(ierr); 930 ierr = PetscFree(qr_work);CHKERRQ(ierr); 931 ierr = PetscFree(gqr_work);CHKERRQ(ierr); 932 ierr = PetscFree(is_primal);CHKERRQ(ierr); 933 ierr = PetscFree(qr_basis);CHKERRQ(ierr); 934 } 935 /* assembling */ 936 ierr = MatAssemblyBegin(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 937 ierr = MatAssemblyEnd(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 938 /* 939 ierr = MatView(pcbddc->ChangeOfBasisMatrix,(PetscViewer)0);CHKERRQ(ierr); 940 */ 941 } 942 /* free workspace no longer needed */ 943 ierr = PetscFree(aux_primal_numbering);CHKERRQ(ierr); 944 ierr = PetscFree(aux_primal_minloc);CHKERRQ(ierr); 945 ierr = PetscFree(temp_indices);CHKERRQ(ierr); 946 ierr = PetscFree(change_basis);CHKERRQ(ierr); 947 ierr = PetscFree(temp_indices_to_constraint);CHKERRQ(ierr); 948 ierr = PetscFree(temp_indices_to_constraint_B);CHKERRQ(ierr); 949 ierr = PetscFree(local_to_B);CHKERRQ(ierr); 950 ierr = PetscFree(temp_quadrature_constraint);CHKERRQ(ierr); 951 PetscFunctionReturn(0); 952 } 953 954 #undef __FUNCT__ 955 #define __FUNCT__ "PCBDDCAnalyzeInterface" 956 PetscErrorCode PCBDDCAnalyzeInterface(PC pc) 957 { 958 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 959 PC_IS *pcis = (PC_IS*)pc->data; 960 Mat_IS *matis = (Mat_IS*)pc->pmat->data; 961 PetscInt bs,ierr,i,vertex_size; 962 PetscViewer viewer=pcbddc->dbg_viewer; 963 964 PetscFunctionBegin; 965 /* Init local Graph struct */ 966 ierr = PCBDDCGraphInit(pcbddc->mat_graph,matis->mapping);CHKERRQ(ierr); 967 968 /* Check validity of the csr graph passed in by the user */ 969 if (pcbddc->mat_graph->nvtxs_csr != pcbddc->mat_graph->nvtxs) { 970 ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr); 971 } 972 /* Set default CSR adjacency of local dofs if not provided by the user with PCBDDCSetLocalAdjacencyGraph */ 973 if (!pcbddc->mat_graph->xadj || !pcbddc->mat_graph->adjncy) { 974 Mat mat_adj; 975 const PetscInt *xadj,*adjncy; 976 PetscBool flg_row=PETSC_TRUE; 977 978 ierr = MatConvert(matis->A,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr); 979 ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr); 980 if (!flg_row) { 981 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__); 982 } 983 ierr = PCBDDCSetLocalAdjacencyGraph(pc,i,xadj,adjncy,PETSC_COPY_VALUES);CHKERRQ(ierr); 984 ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr); 985 if (!flg_row) { 986 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 987 } 988 ierr = MatDestroy(&mat_adj);CHKERRQ(ierr); 989 } 990 991 /* Set default dofs' splitting if no information has been provided by the user with PCBDDCSetDofsSplitting */ 992 vertex_size = 1; 993 if (!pcbddc->n_ISForDofs) { 994 IS *custom_ISForDofs; 995 996 ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 997 ierr = PetscMalloc(bs*sizeof(IS),&custom_ISForDofs);CHKERRQ(ierr); 998 for (i=0;i<bs;i++) { 999 ierr = ISCreateStride(PETSC_COMM_SELF,pcis->n/bs,i,bs,&custom_ISForDofs[i]);CHKERRQ(ierr); 1000 } 1001 ierr = PCBDDCSetDofsSplitting(pc,bs,custom_ISForDofs);CHKERRQ(ierr); 1002 /* remove my references to IS objects */ 1003 for (i=0;i<bs;i++) { 1004 ierr = ISDestroy(&custom_ISForDofs[i]);CHKERRQ(ierr); 1005 } 1006 ierr = PetscFree(custom_ISForDofs);CHKERRQ(ierr); 1007 } else { /* mat block size as vertex size (used for elasticity) */ 1008 ierr = MatGetBlockSize(matis->A,&vertex_size);CHKERRQ(ierr); 1009 } 1010 1011 /* Setup of Graph */ 1012 ierr = PCBDDCGraphSetUp(pcbddc->mat_graph,vertex_size,pcbddc->NeumannBoundaries,pcbddc->DirichletBoundaries,pcbddc->n_ISForDofs,pcbddc->ISForDofs,pcbddc->user_primal_vertices); 1013 1014 /* Graph's connected components analysis */ 1015 ierr = PCBDDCGraphComputeConnectedComponents(pcbddc->mat_graph);CHKERRQ(ierr); 1016 1017 /* print some info to stdout */ 1018 if (pcbddc->dbg_flag) { 1019 ierr = PCBDDCGraphASCIIView(pcbddc->mat_graph,pcbddc->dbg_flag,viewer); 1020 } 1021 PetscFunctionReturn(0); 1022 } 1023 1024 #undef __FUNCT__ 1025 #define __FUNCT__ "PCBDDCGetPrimalVerticesLocalIdx" 1026 PetscErrorCode PCBDDCGetPrimalVerticesLocalIdx(PC pc, PetscInt *n_vertices, PetscInt *vertices_idx[]) 1027 { 1028 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 1029 PetscInt *vertices,*row_cmat_indices,n,i,size_of_constraint,local_primal_size; 1030 PetscErrorCode ierr; 1031 1032 PetscFunctionBegin; 1033 n = 0; 1034 vertices = 0; 1035 if (pcbddc->ConstraintMatrix) { 1036 ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&i);CHKERRQ(ierr); 1037 for (i=0;i<local_primal_size;i++) { 1038 ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 1039 if (size_of_constraint == 1) n++; 1040 ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 1041 } 1042 ierr = PetscMalloc(n*sizeof(PetscInt),&vertices);CHKERRQ(ierr); 1043 n = 0; 1044 for (i=0;i<local_primal_size;i++) { 1045 ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1046 if (size_of_constraint == 1) { 1047 vertices[n++]=row_cmat_indices[0]; 1048 } 1049 ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1050 } 1051 } 1052 *n_vertices = n; 1053 *vertices_idx = vertices; 1054 PetscFunctionReturn(0); 1055 } 1056 1057 #undef __FUNCT__ 1058 #define __FUNCT__ "PCBDDCGetPrimalConstraintsLocalIdx" 1059 PetscErrorCode PCBDDCGetPrimalConstraintsLocalIdx(PC pc, PetscInt *n_constraints, PetscInt *constraints_idx[]) 1060 { 1061 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 1062 PetscInt *constraints_index,*row_cmat_indices,*row_cmat_global_indices; 1063 PetscInt n,i,j,size_of_constraint,local_primal_size,local_size,max_size_of_constraint,min_index,min_loc; 1064 PetscBool *touched; 1065 PetscErrorCode ierr; 1066 1067 PetscFunctionBegin; 1068 n = 0; 1069 constraints_index = 0; 1070 if (pcbddc->ConstraintMatrix) { 1071 ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&local_size);CHKERRQ(ierr); 1072 max_size_of_constraint = 0; 1073 for (i=0;i<local_primal_size;i++) { 1074 ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 1075 if (size_of_constraint > 1) { 1076 n++; 1077 } 1078 max_size_of_constraint = PetscMax(size_of_constraint,max_size_of_constraint); 1079 ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 1080 } 1081 ierr = PetscMalloc(n*sizeof(PetscInt),&constraints_index);CHKERRQ(ierr); 1082 ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&row_cmat_global_indices);CHKERRQ(ierr); 1083 ierr = PetscMalloc(local_size*sizeof(PetscBool),&touched);CHKERRQ(ierr); 1084 ierr = PetscMemzero(touched,local_size*sizeof(PetscBool));CHKERRQ(ierr); 1085 n = 0; 1086 for (i=0;i<local_primal_size;i++) { 1087 ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1088 if (size_of_constraint > 1) { 1089 ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,row_cmat_indices,row_cmat_global_indices);CHKERRQ(ierr); 1090 min_index = row_cmat_global_indices[0]; 1091 min_loc = 0; 1092 for (j=1;j<size_of_constraint;j++) { 1093 /* there can be more than one constraint on a single connected component */ 1094 if (min_index > row_cmat_global_indices[j] && !touched[row_cmat_indices[j]]) { 1095 min_index = row_cmat_global_indices[j]; 1096 min_loc = j; 1097 } 1098 } 1099 touched[row_cmat_indices[min_loc]] = PETSC_TRUE; 1100 constraints_index[n++] = row_cmat_indices[min_loc]; 1101 } 1102 ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1103 } 1104 } 1105 ierr = PetscFree(touched);CHKERRQ(ierr); 1106 ierr = PetscFree(row_cmat_global_indices);CHKERRQ(ierr); 1107 *n_constraints = n; 1108 *constraints_idx = constraints_index; 1109 PetscFunctionReturn(0); 1110 } 1111 1112 /* the next two functions has been adapted from pcis.c */ 1113 #undef __FUNCT__ 1114 #define __FUNCT__ "PCBDDCApplySchur" 1115 PetscErrorCode PCBDDCApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 1116 { 1117 PetscErrorCode ierr; 1118 PC_IS *pcis = (PC_IS*)(pc->data); 1119 1120 PetscFunctionBegin; 1121 if (!vec2_B) { vec2_B = v; } 1122 ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 1123 ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr); 1124 ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 1125 ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr); 1126 ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 1127 PetscFunctionReturn(0); 1128 } 1129 1130 #undef __FUNCT__ 1131 #define __FUNCT__ "PCBDDCApplySchurTranspose" 1132 PetscErrorCode PCBDDCApplySchurTranspose(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 1133 { 1134 PetscErrorCode ierr; 1135 PC_IS *pcis = (PC_IS*)(pc->data); 1136 1137 PetscFunctionBegin; 1138 if (!vec2_B) { vec2_B = v; } 1139 ierr = MatMultTranspose(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 1140 ierr = MatMultTranspose(pcis->A_BI,v,vec1_D);CHKERRQ(ierr); 1141 ierr = KSPSolveTranspose(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 1142 ierr = MatMultTranspose(pcis->A_IB,vec2_D,vec2_B);CHKERRQ(ierr); 1143 ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 1144 PetscFunctionReturn(0); 1145 } 1146 1147 #undef __FUNCT__ 1148 #define __FUNCT__ "PCBDDCSubsetNumbering" 1149 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[]) 1150 { 1151 Vec local_vec,global_vec; 1152 IS seqis,paris; 1153 VecScatter scatter_ctx; 1154 PetscScalar *array; 1155 PetscInt *temp_global_dofs; 1156 PetscScalar globalsum; 1157 PetscInt i,j,s; 1158 PetscInt nlocals,first_index,old_index,max_local; 1159 PetscMPIInt rank_prec_comm,size_prec_comm,max_global; 1160 PetscMPIInt *dof_sizes,*dof_displs; 1161 PetscBool first_found; 1162 PetscErrorCode ierr; 1163 1164 PetscFunctionBegin; 1165 /* mpi buffers */ 1166 MPI_Comm_size(comm,&size_prec_comm); 1167 MPI_Comm_rank(comm,&rank_prec_comm); 1168 j = ( !rank_prec_comm ? size_prec_comm : 0); 1169 ierr = PetscMalloc(j*sizeof(*dof_sizes),&dof_sizes);CHKERRQ(ierr); 1170 ierr = PetscMalloc(j*sizeof(*dof_displs),&dof_displs);CHKERRQ(ierr); 1171 /* get maximum size of subset */ 1172 ierr = PetscMalloc(n_local_dofs*sizeof(PetscInt),&temp_global_dofs);CHKERRQ(ierr); 1173 ierr = ISLocalToGlobalMappingApply(l2gmap,n_local_dofs,local_dofs,temp_global_dofs);CHKERRQ(ierr); 1174 max_local = 0; 1175 if (n_local_dofs) { 1176 max_local = temp_global_dofs[0]; 1177 for (i=1;i<n_local_dofs;i++) { 1178 if (max_local < temp_global_dofs[i] ) { 1179 max_local = temp_global_dofs[i]; 1180 } 1181 } 1182 } 1183 ierr = MPI_Allreduce(&max_local,&max_global,1,MPIU_INT,MPI_MAX,comm); 1184 max_global++; 1185 max_local = 0; 1186 if (n_local_dofs) { 1187 max_local = local_dofs[0]; 1188 for (i=1;i<n_local_dofs;i++) { 1189 if (max_local < local_dofs[i] ) { 1190 max_local = local_dofs[i]; 1191 } 1192 } 1193 } 1194 max_local++; 1195 /* allocate workspace */ 1196 ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr); 1197 ierr = VecSetSizes(local_vec,PETSC_DECIDE,max_local);CHKERRQ(ierr); 1198 ierr = VecSetType(local_vec,VECSEQ);CHKERRQ(ierr); 1199 ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr); 1200 ierr = VecSetSizes(global_vec,PETSC_DECIDE,max_global);CHKERRQ(ierr); 1201 ierr = VecSetType(global_vec,VECMPI);CHKERRQ(ierr); 1202 /* create scatter */ 1203 ierr = ISCreateGeneral(PETSC_COMM_SELF,n_local_dofs,local_dofs,PETSC_COPY_VALUES,&seqis);CHKERRQ(ierr); 1204 ierr = ISCreateGeneral(comm,n_local_dofs,temp_global_dofs,PETSC_COPY_VALUES,&paris);CHKERRQ(ierr); 1205 ierr = VecScatterCreate(local_vec,seqis,global_vec,paris,&scatter_ctx);CHKERRQ(ierr); 1206 ierr = ISDestroy(&seqis);CHKERRQ(ierr); 1207 ierr = ISDestroy(&paris);CHKERRQ(ierr); 1208 /* init array */ 1209 ierr = VecSet(global_vec,0.0);CHKERRQ(ierr); 1210 ierr = VecSet(local_vec,0.0);CHKERRQ(ierr); 1211 ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 1212 if (local_dofs_mult) { 1213 for (i=0;i<n_local_dofs;i++) { 1214 array[local_dofs[i]]=(PetscScalar)local_dofs_mult[i]; 1215 } 1216 } else { 1217 for (i=0;i<n_local_dofs;i++) { 1218 array[local_dofs[i]]=1.0; 1219 } 1220 } 1221 ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 1222 /* scatter into global vec and get total number of global dofs */ 1223 ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1224 ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1225 ierr = VecSum(global_vec,&globalsum);CHKERRQ(ierr); 1226 *n_global_subset = (PetscInt)PetscRealPart(globalsum); 1227 /* Fill global_vec with cumulative function for global numbering */ 1228 ierr = VecGetArray(global_vec,&array);CHKERRQ(ierr); 1229 ierr = VecGetLocalSize(global_vec,&s);CHKERRQ(ierr); 1230 nlocals = 0; 1231 first_index = -1; 1232 first_found = PETSC_FALSE; 1233 for (i=0;i<s;i++) { 1234 if (!first_found && PetscRealPart(array[i]) > 0.0) { 1235 first_found = PETSC_TRUE; 1236 first_index = i; 1237 } 1238 nlocals += (PetscInt)PetscRealPart(array[i]); 1239 } 1240 ierr = MPI_Gather(&nlocals,1,MPIU_INT,dof_sizes,1,MPIU_INT,0,comm);CHKERRQ(ierr); 1241 if (!rank_prec_comm) { 1242 dof_displs[0]=0; 1243 for (i=1;i<size_prec_comm;i++) { 1244 dof_displs[i] = dof_displs[i-1]+dof_sizes[i-1]; 1245 } 1246 } 1247 ierr = MPI_Scatter(dof_displs,1,MPIU_INT,&nlocals,1,MPIU_INT,0,comm);CHKERRQ(ierr); 1248 if (first_found) { 1249 array[first_index] += (PetscScalar)nlocals; 1250 old_index = first_index; 1251 for (i=first_index+1;i<s;i++) { 1252 if (PetscRealPart(array[i]) > 0.0) { 1253 array[i] += array[old_index]; 1254 old_index = i; 1255 } 1256 } 1257 } 1258 ierr = VecRestoreArray(global_vec,&array);CHKERRQ(ierr); 1259 ierr = VecSet(local_vec,0.0);CHKERRQ(ierr); 1260 ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1261 ierr = VecScatterEnd (scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1262 /* get global ordering of local dofs */ 1263 ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 1264 if (local_dofs_mult) { 1265 for (i=0;i<n_local_dofs;i++) { 1266 temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-local_dofs_mult[i]; 1267 } 1268 } else { 1269 for (i=0;i<n_local_dofs;i++) { 1270 temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-1; 1271 } 1272 } 1273 ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 1274 /* free workspace */ 1275 ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr); 1276 ierr = VecDestroy(&local_vec);CHKERRQ(ierr); 1277 ierr = VecDestroy(&global_vec);CHKERRQ(ierr); 1278 ierr = PetscFree(dof_sizes);CHKERRQ(ierr); 1279 ierr = PetscFree(dof_displs);CHKERRQ(ierr); 1280 /* return pointer to global ordering of local dofs */ 1281 *global_numbering_subset = temp_global_dofs; 1282 PetscFunctionReturn(0); 1283 } 1284