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