1674ae819SStefano Zampini #include "bddc.h" 2674ae819SStefano Zampini #include "bddcprivate.h" 3674ae819SStefano Zampini #include <petscblaslapack.h> 4674ae819SStefano Zampini 5674ae819SStefano Zampini #undef __FUNCT__ 6674ae819SStefano Zampini #define __FUNCT__ "PCBDDCResetCustomization" 7674ae819SStefano Zampini PetscErrorCode PCBDDCResetCustomization(PC pc) 8674ae819SStefano Zampini { 9674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 10674ae819SStefano Zampini PetscInt i; 11674ae819SStefano Zampini PetscErrorCode ierr; 12674ae819SStefano Zampini 13674ae819SStefano Zampini PetscFunctionBegin; 14674ae819SStefano Zampini ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr); 15674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 16674ae819SStefano Zampini ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 17674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 18674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 19674ae819SStefano Zampini for (i=0;i<pcbddc->n_ISForDofs;i++) { 20674ae819SStefano Zampini ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 21674ae819SStefano Zampini } 22674ae819SStefano Zampini ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 23674ae819SStefano Zampini PetscFunctionReturn(0); 24674ae819SStefano Zampini } 25674ae819SStefano Zampini 26674ae819SStefano Zampini #undef __FUNCT__ 27674ae819SStefano Zampini #define __FUNCT__ "PCBDDCResetTopography" 28674ae819SStefano Zampini PetscErrorCode PCBDDCResetTopography(PC pc) 29674ae819SStefano Zampini { 30674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 31674ae819SStefano Zampini PetscErrorCode ierr; 32674ae819SStefano Zampini 33674ae819SStefano Zampini PetscFunctionBegin; 34674ae819SStefano Zampini ierr = MatDestroy(&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr); 35674ae819SStefano Zampini ierr = MatDestroy(&pcbddc->ConstraintMatrix);CHKERRQ(ierr); 36674ae819SStefano Zampini ierr = PCBDDCGraphReset(pcbddc->mat_graph);CHKERRQ(ierr); 37674ae819SStefano Zampini PetscFunctionReturn(0); 38674ae819SStefano Zampini } 39674ae819SStefano Zampini 40674ae819SStefano Zampini #undef __FUNCT__ 41674ae819SStefano Zampini #define __FUNCT__ "PCBDDCResetSolvers" 42674ae819SStefano Zampini PetscErrorCode PCBDDCResetSolvers(PC pc) 43674ae819SStefano Zampini { 44674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 45674ae819SStefano Zampini PetscErrorCode ierr; 46674ae819SStefano Zampini 47674ae819SStefano Zampini PetscFunctionBegin; 48674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 49674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 50674ae819SStefano Zampini ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 51674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->coarse_vec);CHKERRQ(ierr); 52674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->coarse_rhs);CHKERRQ(ierr); 53674ae819SStefano Zampini ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr); 54674ae819SStefano Zampini ierr = MatDestroy(&pcbddc->coarse_mat);CHKERRQ(ierr); 55674ae819SStefano Zampini ierr = MatDestroy(&pcbddc->coarse_phi_B);CHKERRQ(ierr); 56674ae819SStefano Zampini ierr = MatDestroy(&pcbddc->coarse_phi_D);CHKERRQ(ierr); 57674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->vec1_P);CHKERRQ(ierr); 58674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->vec1_C);CHKERRQ(ierr); 59674ae819SStefano Zampini ierr = MatDestroy(&pcbddc->local_auxmat1);CHKERRQ(ierr); 60674ae819SStefano Zampini ierr = MatDestroy(&pcbddc->local_auxmat2);CHKERRQ(ierr); 61674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->vec1_R);CHKERRQ(ierr); 62674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->vec2_R);CHKERRQ(ierr); 63674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->vec4_D);CHKERRQ(ierr); 64674ae819SStefano Zampini ierr = VecScatterDestroy(&pcbddc->R_to_B);CHKERRQ(ierr); 65674ae819SStefano Zampini ierr = VecScatterDestroy(&pcbddc->R_to_D);CHKERRQ(ierr); 66674ae819SStefano Zampini ierr = VecScatterDestroy(&pcbddc->coarse_loc_to_glob);CHKERRQ(ierr); 67674ae819SStefano Zampini ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr); 68674ae819SStefano Zampini ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr); 69674ae819SStefano Zampini ierr = PetscFree(pcbddc->local_primal_indices);CHKERRQ(ierr); 70674ae819SStefano Zampini ierr = PetscFree(pcbddc->replicated_local_primal_indices);CHKERRQ(ierr); 71674ae819SStefano Zampini ierr = PetscFree(pcbddc->replicated_local_primal_values);CHKERRQ(ierr); 72674ae819SStefano Zampini ierr = PetscFree(pcbddc->local_primal_displacements);CHKERRQ(ierr); 73674ae819SStefano Zampini ierr = PetscFree(pcbddc->local_primal_sizes);CHKERRQ(ierr); 74674ae819SStefano Zampini PetscFunctionReturn(0); 75674ae819SStefano Zampini } 76674ae819SStefano Zampini 77674ae819SStefano Zampini #undef __FUNCT__ 78674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSolveSaddlePoint" 79674ae819SStefano Zampini static PetscErrorCode PCBDDCSolveSaddlePoint(PC pc) 80674ae819SStefano Zampini { 81674ae819SStefano Zampini PetscErrorCode ierr; 82674ae819SStefano Zampini PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 83674ae819SStefano Zampini 84674ae819SStefano Zampini PetscFunctionBegin; 85674ae819SStefano Zampini ierr = KSPSolve(pcbddc->ksp_R,pcbddc->vec1_R,pcbddc->vec2_R);CHKERRQ(ierr); 86674ae819SStefano Zampini if (pcbddc->local_auxmat1) { 87674ae819SStefano Zampini ierr = MatMult(pcbddc->local_auxmat1,pcbddc->vec2_R,pcbddc->vec1_C);CHKERRQ(ierr); 88674ae819SStefano Zampini ierr = MatMultAdd(pcbddc->local_auxmat2,pcbddc->vec1_C,pcbddc->vec2_R,pcbddc->vec2_R);CHKERRQ(ierr); 89674ae819SStefano Zampini } 90674ae819SStefano Zampini PetscFunctionReturn(0); 91674ae819SStefano Zampini } 92674ae819SStefano Zampini 93674ae819SStefano Zampini #undef __FUNCT__ 94674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplyInterfacePreconditioner" 95674ae819SStefano Zampini PetscErrorCode PCBDDCApplyInterfacePreconditioner(PC pc) 96674ae819SStefano Zampini { 97674ae819SStefano Zampini PetscErrorCode ierr; 98674ae819SStefano Zampini PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 99674ae819SStefano Zampini PC_IS* pcis = (PC_IS*) (pc->data); 100674ae819SStefano Zampini const PetscScalar zero = 0.0; 101674ae819SStefano Zampini 102674ae819SStefano Zampini PetscFunctionBegin; 103674ae819SStefano Zampini /* Application of PHI^T */ 104674ae819SStefano Zampini ierr = MatMultTranspose(pcbddc->coarse_phi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr); 105674ae819SStefano Zampini if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_phi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); } 106674ae819SStefano Zampini 107674ae819SStefano Zampini /* Scatter data of coarse_rhs */ 108674ae819SStefano Zampini if (pcbddc->coarse_rhs) { ierr = VecSet(pcbddc->coarse_rhs,zero);CHKERRQ(ierr); } 109674ae819SStefano Zampini ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 110674ae819SStefano Zampini 111674ae819SStefano Zampini /* Local solution on R nodes */ 112674ae819SStefano Zampini ierr = VecSet(pcbddc->vec1_R,zero);CHKERRQ(ierr); 113674ae819SStefano Zampini ierr = VecScatterBegin(pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 114674ae819SStefano Zampini ierr = VecScatterEnd (pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 115674ae819SStefano Zampini if (pcbddc->inexact_prec_type) { 116674ae819SStefano Zampini ierr = VecScatterBegin(pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 117674ae819SStefano Zampini ierr = VecScatterEnd (pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 118674ae819SStefano Zampini } 119674ae819SStefano Zampini ierr = PCBDDCSolveSaddlePoint(pc);CHKERRQ(ierr); 120674ae819SStefano Zampini ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr); 121674ae819SStefano Zampini ierr = VecScatterBegin(pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 122674ae819SStefano Zampini ierr = VecScatterEnd (pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 123674ae819SStefano Zampini if (pcbddc->inexact_prec_type) { 124674ae819SStefano Zampini ierr = VecScatterBegin(pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 125674ae819SStefano Zampini ierr = VecScatterEnd (pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 126674ae819SStefano Zampini } 127674ae819SStefano Zampini 128674ae819SStefano Zampini /* Coarse solution */ 129674ae819SStefano Zampini ierr = PCBDDCScatterCoarseDataEnd(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 130674ae819SStefano Zampini if (pcbddc->coarse_rhs) { /* TODO remove null space when doing multilevel */ 131674ae819SStefano Zampini ierr = KSPSolve(pcbddc->coarse_ksp,pcbddc->coarse_rhs,pcbddc->coarse_vec);CHKERRQ(ierr); 132674ae819SStefano Zampini } 133674ae819SStefano Zampini ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 134674ae819SStefano Zampini ierr = PCBDDCScatterCoarseDataEnd (pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 135674ae819SStefano Zampini 136674ae819SStefano Zampini /* Sum contributions from two levels */ 137674ae819SStefano Zampini ierr = MatMultAdd(pcbddc->coarse_phi_B,pcbddc->vec1_P,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr); 138674ae819SStefano Zampini if (pcbddc->inexact_prec_type) { ierr = MatMultAdd(pcbddc->coarse_phi_D,pcbddc->vec1_P,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 139674ae819SStefano Zampini PetscFunctionReturn(0); 140674ae819SStefano Zampini } 141674ae819SStefano Zampini 142674ae819SStefano Zampini #undef __FUNCT__ 143674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScatterCoarseDataBegin" 144674ae819SStefano Zampini PetscErrorCode PCBDDCScatterCoarseDataBegin(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode) 145674ae819SStefano Zampini { 146674ae819SStefano Zampini PetscErrorCode ierr; 147674ae819SStefano Zampini PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 148674ae819SStefano Zampini 149674ae819SStefano Zampini PetscFunctionBegin; 150674ae819SStefano Zampini switch (pcbddc->coarse_communications_type) { 151674ae819SStefano Zampini case SCATTERS_BDDC: 152674ae819SStefano Zampini ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr); 153674ae819SStefano Zampini break; 154674ae819SStefano Zampini case GATHERS_BDDC: 155674ae819SStefano Zampini break; 156674ae819SStefano Zampini } 157674ae819SStefano Zampini PetscFunctionReturn(0); 158674ae819SStefano Zampini } 159674ae819SStefano Zampini 160674ae819SStefano Zampini #undef __FUNCT__ 161674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScatterCoarseDataEnd" 162674ae819SStefano Zampini PetscErrorCode PCBDDCScatterCoarseDataEnd(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode) 163674ae819SStefano Zampini { 164674ae819SStefano Zampini PetscErrorCode ierr; 165674ae819SStefano Zampini PC_BDDC* pcbddc = (PC_BDDC*)(pc->data); 166674ae819SStefano Zampini PetscScalar* array_to; 167674ae819SStefano Zampini PetscScalar* array_from; 168674ae819SStefano Zampini MPI_Comm comm; 169674ae819SStefano Zampini PetscInt i; 170674ae819SStefano Zampini 171674ae819SStefano Zampini PetscFunctionBegin; 172674ae819SStefano Zampini ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 173674ae819SStefano Zampini switch (pcbddc->coarse_communications_type) { 174674ae819SStefano Zampini case SCATTERS_BDDC: 175674ae819SStefano Zampini ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr); 176674ae819SStefano Zampini break; 177674ae819SStefano Zampini case GATHERS_BDDC: 178674ae819SStefano Zampini if (vec_from) { 179674ae819SStefano Zampini ierr = VecGetArray(vec_from,&array_from);CHKERRQ(ierr); 180674ae819SStefano Zampini } 181674ae819SStefano Zampini if (vec_to) { 182674ae819SStefano Zampini ierr = VecGetArray(vec_to,&array_to);CHKERRQ(ierr); 183674ae819SStefano Zampini } 184674ae819SStefano Zampini switch(pcbddc->coarse_problem_type){ 185674ae819SStefano Zampini case SEQUENTIAL_BDDC: 186674ae819SStefano Zampini if (smode == SCATTER_FORWARD) { 187674ae819SStefano Zampini 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); 188674ae819SStefano Zampini if (vec_to) { 189674ae819SStefano Zampini if (imode == ADD_VALUES) { 190674ae819SStefano Zampini for (i=0;i<pcbddc->replicated_primal_size;i++) { 191674ae819SStefano Zampini array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i]; 192674ae819SStefano Zampini } 193674ae819SStefano Zampini } else { 194674ae819SStefano Zampini for (i=0;i<pcbddc->replicated_primal_size;i++) { 195674ae819SStefano Zampini array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i]; 196674ae819SStefano Zampini } 197674ae819SStefano Zampini } 198674ae819SStefano Zampini } 199674ae819SStefano Zampini } else { 200674ae819SStefano Zampini if (vec_from) { 201674ae819SStefano Zampini if (imode == ADD_VALUES) { 202674ae819SStefano Zampini MPI_Comm vec_from_comm; 203674ae819SStefano Zampini ierr = PetscObjectGetComm((PetscObject)(vec_from),&vec_from_comm);CHKERRQ(ierr); 204674ae819SStefano Zampini 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); 205674ae819SStefano Zampini } 206674ae819SStefano Zampini for (i=0;i<pcbddc->replicated_primal_size;i++) { 207674ae819SStefano Zampini pcbddc->replicated_local_primal_values[i]=array_from[pcbddc->replicated_local_primal_indices[i]]; 208674ae819SStefano Zampini } 209674ae819SStefano Zampini } 210674ae819SStefano Zampini 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); 211674ae819SStefano Zampini } 212674ae819SStefano Zampini break; 213674ae819SStefano Zampini case REPLICATED_BDDC: 214674ae819SStefano Zampini if (smode == SCATTER_FORWARD) { 215674ae819SStefano Zampini 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); 216674ae819SStefano Zampini if (imode == ADD_VALUES) { 217674ae819SStefano Zampini for (i=0;i<pcbddc->replicated_primal_size;i++) { 218674ae819SStefano Zampini array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i]; 219674ae819SStefano Zampini } 220674ae819SStefano Zampini } else { 221674ae819SStefano Zampini for (i=0;i<pcbddc->replicated_primal_size;i++) { 222674ae819SStefano Zampini array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i]; 223674ae819SStefano Zampini } 224674ae819SStefano Zampini } 225674ae819SStefano Zampini } else { /* no communications needed for SCATTER_REVERSE since needed data is already present */ 226674ae819SStefano Zampini if (imode == ADD_VALUES) { 227674ae819SStefano Zampini for (i=0;i<pcbddc->local_primal_size;i++) { 228674ae819SStefano Zampini array_to[i]+=array_from[pcbddc->local_primal_indices[i]]; 229674ae819SStefano Zampini } 230674ae819SStefano Zampini } else { 231674ae819SStefano Zampini for (i=0;i<pcbddc->local_primal_size;i++) { 232674ae819SStefano Zampini array_to[i]=array_from[pcbddc->local_primal_indices[i]]; 233674ae819SStefano Zampini } 234674ae819SStefano Zampini } 235674ae819SStefano Zampini } 236674ae819SStefano Zampini break; 237674ae819SStefano Zampini case MULTILEVEL_BDDC: 238674ae819SStefano Zampini break; 239674ae819SStefano Zampini case PARALLEL_BDDC: 240674ae819SStefano Zampini break; 241674ae819SStefano Zampini } 242674ae819SStefano Zampini if (vec_from) { 243674ae819SStefano Zampini ierr = VecRestoreArray(vec_from,&array_from);CHKERRQ(ierr); 244674ae819SStefano Zampini } 245674ae819SStefano Zampini if (vec_to) { 246674ae819SStefano Zampini ierr = VecRestoreArray(vec_to,&array_to);CHKERRQ(ierr); 247674ae819SStefano Zampini } 248674ae819SStefano Zampini break; 249674ae819SStefano Zampini } 250674ae819SStefano Zampini PetscFunctionReturn(0); 251674ae819SStefano Zampini } 252674ae819SStefano Zampini 253674ae819SStefano Zampini #undef __FUNCT__ 254674ae819SStefano Zampini #define __FUNCT__ "PCBDDCConstraintsSetUp" 255674ae819SStefano Zampini PetscErrorCode PCBDDCConstraintsSetUp(PC pc) 256674ae819SStefano Zampini { 257674ae819SStefano Zampini PetscErrorCode ierr; 258674ae819SStefano Zampini PC_IS* pcis = (PC_IS*)(pc->data); 259674ae819SStefano Zampini PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 260674ae819SStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 261674ae819SStefano Zampini PetscInt *nnz,*is_indices; 262674ae819SStefano Zampini PetscScalar *temp_quadrature_constraint; 263674ae819SStefano Zampini PetscInt *temp_indices,*temp_indices_to_constraint,*temp_indices_to_constraint_B,*local_to_B; 264674ae819SStefano Zampini PetscInt local_primal_size,i,j,k,total_counts,max_size_of_constraint; 265674ae819SStefano Zampini PetscInt n_vertices,size_of_constraint; 266674ae819SStefano Zampini PetscScalar quad_value; 267674ae819SStefano Zampini PetscBool nnsp_has_cnst=PETSC_FALSE,use_nnsp_true=pcbddc->use_nnsp_true; 268674ae819SStefano Zampini PetscInt nnsp_size=0,nnsp_addone=0,temp_constraints,temp_start_ptr,n_ISForFaces,n_ISForEdges; 269674ae819SStefano Zampini IS *used_IS,ISForVertices,*ISForFaces,*ISForEdges; 270674ae819SStefano Zampini MatType impMatType=MATSEQAIJ; 271674ae819SStefano Zampini PetscBLASInt Bs,Bt,lwork,lierr; 272674ae819SStefano Zampini PetscReal tol=1.0e-8; 273674ae819SStefano Zampini MatNullSpace nearnullsp; 274674ae819SStefano Zampini const Vec *nearnullvecs; 275674ae819SStefano Zampini Vec *localnearnullsp; 276674ae819SStefano Zampini PetscScalar *work,*temp_basis,*array_vector,*correlation_mat; 277674ae819SStefano Zampini PetscReal *rwork,*singular_vals; 278674ae819SStefano Zampini PetscBLASInt Bone=1,*ipiv; 279674ae819SStefano Zampini Vec temp_vec; 280674ae819SStefano Zampini Mat temp_mat; 281674ae819SStefano Zampini KSP temp_ksp; 282674ae819SStefano Zampini PC temp_pc; 283674ae819SStefano Zampini PetscInt s,start_constraint,dual_dofs; 284674ae819SStefano Zampini PetscBool compute_submatrix,useksp=PETSC_FALSE; 285674ae819SStefano Zampini PetscInt *aux_primal_permutation,*aux_primal_numbering; 286674ae819SStefano Zampini PetscBool boolforface,*change_basis; 287674ae819SStefano Zampini /* some ugly conditional declarations */ 288674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD) 289674ae819SStefano Zampini PetscScalar one=1.0,zero=0.0; 290674ae819SStefano Zampini PetscInt ii; 291674ae819SStefano Zampini PetscScalar *singular_vectors; 292674ae819SStefano Zampini PetscBLASInt *iwork,*ifail; 293674ae819SStefano Zampini PetscReal dummy_real,abs_tol; 294674ae819SStefano Zampini PetscBLASInt eigs_found; 295674ae819SStefano Zampini #endif 296674ae819SStefano Zampini PetscBLASInt dummy_int; 297674ae819SStefano Zampini PetscScalar dummy_scalar; 298674ae819SStefano Zampini PetscBool used_vertex,get_faces,get_edges,get_vertices; 299674ae819SStefano Zampini 300674ae819SStefano Zampini PetscFunctionBegin; 301674ae819SStefano Zampini /* Get index sets for faces, edges and vertices from graph */ 302674ae819SStefano Zampini get_faces = PETSC_TRUE; 303674ae819SStefano Zampini get_edges = PETSC_TRUE; 304674ae819SStefano Zampini get_vertices = PETSC_TRUE; 305674ae819SStefano Zampini if (pcbddc->vertices_flag) { 306674ae819SStefano Zampini get_faces = PETSC_FALSE; 307674ae819SStefano Zampini get_edges = PETSC_FALSE; 308674ae819SStefano Zampini } 309674ae819SStefano Zampini if (pcbddc->constraints_flag) { 310674ae819SStefano Zampini get_vertices = PETSC_FALSE; 311674ae819SStefano Zampini } 312674ae819SStefano Zampini if (pcbddc->faces_flag) { 313674ae819SStefano Zampini get_edges = PETSC_FALSE; 314674ae819SStefano Zampini } 315674ae819SStefano Zampini if (pcbddc->edges_flag) { 316674ae819SStefano Zampini get_faces = PETSC_FALSE; 317674ae819SStefano Zampini } 318674ae819SStefano Zampini /* default */ 319674ae819SStefano Zampini if (!get_faces && !get_edges && !get_vertices) { 320674ae819SStefano Zampini get_vertices = PETSC_TRUE; 321674ae819SStefano Zampini } 322674ae819SStefano Zampini ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,get_faces,get_edges,get_vertices,&n_ISForFaces,&ISForFaces,&n_ISForEdges,&ISForEdges,&ISForVertices); 323674ae819SStefano Zampini if (pcbddc->dbg_flag) { 324674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr); 325674ae819SStefano Zampini i = 0; 326674ae819SStefano Zampini if (ISForVertices) { 327674ae819SStefano Zampini ierr = ISGetSize(ISForVertices,&i);CHKERRQ(ierr); 328674ae819SStefano Zampini } 329674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate vertices\n",PetscGlobalRank,i);CHKERRQ(ierr); 330674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate faces\n",PetscGlobalRank,n_ISForFaces);CHKERRQ(ierr); 331674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate edges\n",PetscGlobalRank,n_ISForEdges);CHKERRQ(ierr); 332674ae819SStefano Zampini ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 333674ae819SStefano Zampini } 334674ae819SStefano Zampini /* check if near null space is attached to global mat */ 335674ae819SStefano Zampini ierr = MatGetNearNullSpace(pc->pmat,&nearnullsp);CHKERRQ(ierr); 336674ae819SStefano Zampini if (nearnullsp) { 337674ae819SStefano Zampini ierr = MatNullSpaceGetVecs(nearnullsp,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr); 338674ae819SStefano Zampini } else { /* if near null space is not provided it uses constants */ 339674ae819SStefano Zampini nnsp_has_cnst = PETSC_TRUE; 340674ae819SStefano Zampini use_nnsp_true = PETSC_TRUE; 341674ae819SStefano Zampini } 342674ae819SStefano Zampini if (nnsp_has_cnst) { 343674ae819SStefano Zampini nnsp_addone = 1; 344674ae819SStefano Zampini } 345674ae819SStefano Zampini /* 346674ae819SStefano Zampini Evaluate maximum storage size needed by the procedure 347674ae819SStefano Zampini - temp_indices will contain start index of each constraint stored as follows 348674ae819SStefano Zampini - temp_indices_to_constraint [temp_indices[i],...,temp[indices[i+1]-1] will contain the indices (in local numbering) on which the constraint acts 349674ae819SStefano Zampini - 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 350674ae819SStefano Zampini - temp_quadrature_constraint [temp_indices[i],...,temp[indices[i+1]-1] will contain the scalars representing the constraint itself 351674ae819SStefano Zampini */ 352674ae819SStefano Zampini total_counts = n_ISForFaces+n_ISForEdges; 353674ae819SStefano Zampini total_counts *= (nnsp_addone+nnsp_size); 354674ae819SStefano Zampini n_vertices = 0; 355674ae819SStefano Zampini if (ISForVertices) { 356674ae819SStefano Zampini ierr = ISGetSize(ISForVertices,&n_vertices);CHKERRQ(ierr); 357674ae819SStefano Zampini } 358674ae819SStefano Zampini total_counts += n_vertices; 359674ae819SStefano Zampini ierr = PetscMalloc((total_counts+1)*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr); 360674ae819SStefano Zampini ierr = PetscMalloc((total_counts+1)*sizeof(PetscBool),&change_basis);CHKERRQ(ierr); 361674ae819SStefano Zampini total_counts = 0; 362674ae819SStefano Zampini max_size_of_constraint = 0; 363674ae819SStefano Zampini for (i=0;i<n_ISForEdges+n_ISForFaces;i++) { 364674ae819SStefano Zampini if (i<n_ISForEdges) { 365674ae819SStefano Zampini used_IS = &ISForEdges[i]; 366674ae819SStefano Zampini } else { 367674ae819SStefano Zampini used_IS = &ISForFaces[i-n_ISForEdges]; 368674ae819SStefano Zampini } 369674ae819SStefano Zampini ierr = ISGetSize(*used_IS,&j);CHKERRQ(ierr); 370674ae819SStefano Zampini total_counts += j; 371674ae819SStefano Zampini max_size_of_constraint = PetscMax(j,max_size_of_constraint); 372674ae819SStefano Zampini } 373674ae819SStefano Zampini total_counts *= (nnsp_addone+nnsp_size); 374674ae819SStefano Zampini total_counts += n_vertices; 375674ae819SStefano Zampini ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&temp_quadrature_constraint);CHKERRQ(ierr); 376674ae819SStefano Zampini ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint);CHKERRQ(ierr); 377674ae819SStefano Zampini ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint_B);CHKERRQ(ierr); 378674ae819SStefano Zampini ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&local_to_B);CHKERRQ(ierr); 379674ae819SStefano Zampini ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 380674ae819SStefano Zampini for (i=0;i<pcis->n;i++) { 381674ae819SStefano Zampini local_to_B[i]=-1; 382674ae819SStefano Zampini } 383674ae819SStefano Zampini for (i=0;i<pcis->n_B;i++) { 384674ae819SStefano Zampini local_to_B[is_indices[i]]=i; 385674ae819SStefano Zampini } 386674ae819SStefano Zampini ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr); 387674ae819SStefano Zampini 388674ae819SStefano Zampini /* First we issue queries to allocate optimal workspace for LAPACKgesvd or LAPACKsyev/LAPACKheev */ 389674ae819SStefano Zampini rwork = 0; 390674ae819SStefano Zampini work = 0; 391674ae819SStefano Zampini singular_vals = 0; 392674ae819SStefano Zampini temp_basis = 0; 393674ae819SStefano Zampini correlation_mat = 0; 394674ae819SStefano Zampini if (!pcbddc->use_nnsp_true) { 395674ae819SStefano Zampini PetscScalar temp_work; 396674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD) 397674ae819SStefano Zampini /* POD */ 398674ae819SStefano Zampini PetscInt max_n; 399674ae819SStefano Zampini max_n = nnsp_addone+nnsp_size; 400674ae819SStefano Zampini /* using some techniques borrowed from Proper Orthogonal Decomposition */ 401674ae819SStefano Zampini ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&correlation_mat);CHKERRQ(ierr); 402674ae819SStefano Zampini ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&singular_vectors);CHKERRQ(ierr); 403674ae819SStefano Zampini ierr = PetscMalloc(max_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr); 404674ae819SStefano Zampini ierr = PetscMalloc(max_size_of_constraint*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr); 405674ae819SStefano Zampini #if defined(PETSC_USE_COMPLEX) 406674ae819SStefano Zampini ierr = PetscMalloc(3*max_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr); 407674ae819SStefano Zampini #endif 408674ae819SStefano Zampini ierr = PetscMalloc(5*max_n*sizeof(PetscBLASInt),&iwork);CHKERRQ(ierr); 409674ae819SStefano Zampini ierr = PetscMalloc(max_n*sizeof(PetscBLASInt),&ifail);CHKERRQ(ierr); 410674ae819SStefano Zampini /* now we evaluate the optimal workspace using query with lwork=-1 */ 411674ae819SStefano Zampini ierr = PetscBLASIntCast(max_n,&Bt);CHKERRQ(ierr); 412674ae819SStefano Zampini lwork=-1; 413674ae819SStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 414674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX) 415674ae819SStefano Zampini abs_tol=1.e-8; 416674ae819SStefano Zampini /* LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,&lierr); */ 417674ae819SStefano Zampini LAPACKsyevx_("V","A","U",&Bt,correlation_mat,&Bt,&dummy_real,&dummy_real,&dummy_int,&dummy_int, 418674ae819SStefano Zampini &abs_tol,&eigs_found,singular_vals,singular_vectors,&Bt,&temp_work,&lwork,iwork,ifail,&lierr); 419674ae819SStefano Zampini #else 420674ae819SStefano Zampini /* LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,rwork,&lierr); */ 421674ae819SStefano Zampini /* LAPACK call is missing here! TODO */ 422674ae819SStefano Zampini SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1"); 423674ae819SStefano Zampini #endif 424674ae819SStefano Zampini if ( lierr ) { 425674ae819SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYEVX Lapack routine %d",(int)lierr); 426674ae819SStefano Zampini } 427674ae819SStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 428674ae819SStefano Zampini #else /* on missing GESVD */ 429674ae819SStefano Zampini /* SVD */ 430674ae819SStefano Zampini PetscInt max_n,min_n; 431674ae819SStefano Zampini max_n = max_size_of_constraint; 432674ae819SStefano Zampini min_n = nnsp_addone+nnsp_size; 433674ae819SStefano Zampini if (max_size_of_constraint < ( nnsp_addone+nnsp_size ) ) { 434674ae819SStefano Zampini min_n = max_size_of_constraint; 435674ae819SStefano Zampini max_n = nnsp_addone+nnsp_size; 436674ae819SStefano Zampini } 437674ae819SStefano Zampini ierr = PetscMalloc(min_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr); 438674ae819SStefano Zampini #if defined(PETSC_USE_COMPLEX) 439674ae819SStefano Zampini ierr = PetscMalloc(5*min_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr); 440674ae819SStefano Zampini #endif 441674ae819SStefano Zampini /* now we evaluate the optimal workspace using query with lwork=-1 */ 442674ae819SStefano Zampini lwork=-1; 443674ae819SStefano Zampini ierr = PetscBLASIntCast(max_n,&Bs);CHKERRQ(ierr); 444674ae819SStefano Zampini ierr = PetscBLASIntCast(min_n,&Bt);CHKERRQ(ierr); 445674ae819SStefano Zampini dummy_int = Bs; 446674ae819SStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 447674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX) 448674ae819SStefano Zampini LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Bs,singular_vals, 449674ae819SStefano Zampini &dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,&temp_work,&lwork,&lierr); 450674ae819SStefano Zampini #else 451674ae819SStefano Zampini LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Bs,singular_vals, 452674ae819SStefano Zampini &dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,&temp_work,&lwork,rwork,&lierr); 453674ae819SStefano Zampini #endif 454674ae819SStefano Zampini if ( lierr ) { 455674ae819SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SVD Lapack routine %d",(int)lierr); 456674ae819SStefano Zampini } 457674ae819SStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 458674ae819SStefano Zampini #endif 459674ae819SStefano Zampini /* Allocate optimal workspace */ 460674ae819SStefano Zampini ierr = PetscBLASIntCast((PetscInt)PetscRealPart(temp_work),&lwork);CHKERRQ(ierr); 461674ae819SStefano Zampini total_counts = (PetscInt)lwork; 462674ae819SStefano Zampini ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&work);CHKERRQ(ierr); 463674ae819SStefano Zampini } 464674ae819SStefano Zampini /* get local part of global near null space vectors */ 465674ae819SStefano Zampini ierr = PetscMalloc(nnsp_size*sizeof(Vec),&localnearnullsp);CHKERRQ(ierr); 466674ae819SStefano Zampini for (k=0;k<nnsp_size;k++) { 467674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_N,&localnearnullsp[k]);CHKERRQ(ierr); 468674ae819SStefano Zampini ierr = VecScatterBegin(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 469674ae819SStefano Zampini ierr = VecScatterEnd(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 470674ae819SStefano Zampini } 471674ae819SStefano Zampini /* Now we can loop on constraining sets */ 472674ae819SStefano Zampini total_counts = 0; 473674ae819SStefano Zampini temp_indices[0] = 0; 474674ae819SStefano Zampini /* vertices */ 475674ae819SStefano Zampini if (ISForVertices) { 476674ae819SStefano Zampini ierr = ISGetIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); 477674ae819SStefano Zampini if (nnsp_has_cnst) { /* consider all vertices */ 478674ae819SStefano Zampini for (i=0;i<n_vertices;i++) { 479674ae819SStefano Zampini temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i]; 480674ae819SStefano Zampini temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]]; 481674ae819SStefano Zampini temp_quadrature_constraint[temp_indices[total_counts]]=1.0; 482674ae819SStefano Zampini temp_indices[total_counts+1]=temp_indices[total_counts]+1; 483674ae819SStefano Zampini change_basis[total_counts]=PETSC_FALSE; 484674ae819SStefano Zampini total_counts++; 485674ae819SStefano Zampini } 486674ae819SStefano Zampini } else { /* consider vertices for which exist at least a localnearnullsp which is not null there */ 487674ae819SStefano Zampini for (i=0;i<n_vertices;i++) { 488674ae819SStefano Zampini used_vertex=PETSC_FALSE; 489674ae819SStefano Zampini k=0; 490674ae819SStefano Zampini while (!used_vertex && k<nnsp_size) { 491674ae819SStefano Zampini ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 492674ae819SStefano Zampini if (PetscAbsScalar(array_vector[is_indices[i]])>0.0) { 493674ae819SStefano Zampini temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i]; 494674ae819SStefano Zampini temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]]; 495674ae819SStefano Zampini temp_quadrature_constraint[temp_indices[total_counts]]=1.0; 496674ae819SStefano Zampini temp_indices[total_counts+1]=temp_indices[total_counts]+1; 497674ae819SStefano Zampini change_basis[total_counts]=PETSC_FALSE; 498674ae819SStefano Zampini total_counts++; 499674ae819SStefano Zampini used_vertex=PETSC_TRUE; 500674ae819SStefano Zampini } 501674ae819SStefano Zampini ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 502674ae819SStefano Zampini k++; 503674ae819SStefano Zampini } 504674ae819SStefano Zampini } 505674ae819SStefano Zampini } 506674ae819SStefano Zampini ierr = ISRestoreIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); 507674ae819SStefano Zampini n_vertices = total_counts; 508674ae819SStefano Zampini } 509674ae819SStefano Zampini /* edges and faces */ 510674ae819SStefano Zampini for (i=0;i<n_ISForEdges+n_ISForFaces;i++) { 511674ae819SStefano Zampini if (i<n_ISForEdges) { 512674ae819SStefano Zampini used_IS = &ISForEdges[i]; 513674ae819SStefano Zampini boolforface = pcbddc->use_change_of_basis; 514674ae819SStefano Zampini } else { 515674ae819SStefano Zampini used_IS = &ISForFaces[i-n_ISForEdges]; 516674ae819SStefano Zampini boolforface = pcbddc->use_change_on_faces; 517674ae819SStefano Zampini } 518674ae819SStefano Zampini temp_constraints = 0; /* zero the number of constraints I have on this conn comp */ 519674ae819SStefano Zampini temp_start_ptr = total_counts; /* need to know the starting index of constraints stored */ 520674ae819SStefano Zampini ierr = ISGetSize(*used_IS,&size_of_constraint);CHKERRQ(ierr); 521674ae819SStefano Zampini ierr = ISGetIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 522674ae819SStefano Zampini if (nnsp_has_cnst) { 523674ae819SStefano Zampini temp_constraints++; 524674ae819SStefano Zampini quad_value = (PetscScalar)(1.0/PetscSqrtReal((PetscReal)size_of_constraint)); 525674ae819SStefano Zampini for (j=0;j<size_of_constraint;j++) { 526674ae819SStefano Zampini temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j]; 527674ae819SStefano Zampini temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]]; 528674ae819SStefano Zampini temp_quadrature_constraint[temp_indices[total_counts]+j]=quad_value; 529674ae819SStefano Zampini } 530674ae819SStefano Zampini temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint; /* store new starting point */ 531674ae819SStefano Zampini change_basis[total_counts]=boolforface; 532674ae819SStefano Zampini total_counts++; 533674ae819SStefano Zampini } 534674ae819SStefano Zampini for (k=0;k<nnsp_size;k++) { 535674ae819SStefano Zampini ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 536674ae819SStefano Zampini for (j=0;j<size_of_constraint;j++) { 537674ae819SStefano Zampini temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j]; 538674ae819SStefano Zampini temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]]; 539674ae819SStefano Zampini temp_quadrature_constraint[temp_indices[total_counts]+j]=array_vector[is_indices[j]]; 540674ae819SStefano Zampini } 541674ae819SStefano Zampini ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr); 542674ae819SStefano Zampini quad_value = 1.0; 543674ae819SStefano Zampini if (use_nnsp_true) { /* check if array is null on the connected component in case use_nnsp_true has been requested */ 544674ae819SStefano Zampini ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr); 545674ae819SStefano Zampini quad_value = BLASasum_(&Bs,&temp_quadrature_constraint[temp_indices[total_counts]],&Bone); 546674ae819SStefano Zampini } 547674ae819SStefano Zampini if (quad_value > 0.0) { /* keep indices and values */ 548674ae819SStefano Zampini temp_constraints++; 549674ae819SStefano Zampini temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint; /* store new starting point */ 550674ae819SStefano Zampini change_basis[total_counts]=boolforface; 551674ae819SStefano Zampini total_counts++; 552674ae819SStefano Zampini } 553674ae819SStefano Zampini } 554674ae819SStefano Zampini ierr = ISRestoreIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 555674ae819SStefano Zampini /* perform SVD on the constraint if use_nnsp_true has not be requested by the user */ 556674ae819SStefano Zampini if (!use_nnsp_true) { 557674ae819SStefano Zampini ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr); 558674ae819SStefano Zampini ierr = PetscBLASIntCast(temp_constraints,&Bt);CHKERRQ(ierr); 559674ae819SStefano Zampini 560674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD) 561674ae819SStefano Zampini ierr = PetscMemzero(correlation_mat,Bt*Bt*sizeof(PetscScalar));CHKERRQ(ierr); 562674ae819SStefano Zampini /* Store upper triangular part of correlation matrix */ 563674ae819SStefano Zampini for (j=0;j<temp_constraints;j++) { 564674ae819SStefano Zampini for (k=0;k<j+1;k++) { 565674ae819SStefano Zampini correlation_mat[j*temp_constraints+k]=BLASdot_(&Bs,&temp_quadrature_constraint[temp_indices[temp_start_ptr+j]],&Bone, 566674ae819SStefano Zampini &temp_quadrature_constraint[temp_indices[temp_start_ptr+k]],&Bone); 567674ae819SStefano Zampini 568674ae819SStefano Zampini } 569674ae819SStefano Zampini } 570674ae819SStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 571674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX) 572674ae819SStefano Zampini /* LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,work,&lwork,&lierr); */ 573674ae819SStefano Zampini LAPACKsyevx_("V","A","U",&Bt,correlation_mat,&Bt,&dummy_real,&dummy_real,&dummy_int,&dummy_int, 574674ae819SStefano Zampini &abs_tol,&eigs_found,singular_vals,singular_vectors,&Bt,work,&lwork,iwork,ifail,&lierr); 575674ae819SStefano Zampini #else 576674ae819SStefano Zampini /* LAPACK call is missing here! TODO */ 577674ae819SStefano Zampini SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1"); 578674ae819SStefano Zampini #endif 579674ae819SStefano Zampini if (lierr) { 580674ae819SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYEVX Lapack routine %d",(int)lierr); 581674ae819SStefano Zampini } 582674ae819SStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 583674ae819SStefano Zampini /* retain eigenvalues greater than tol: note that lapack SYEV gives eigs in ascending order */ 584674ae819SStefano Zampini j=0; 585674ae819SStefano Zampini while (j < Bt && singular_vals[j] < tol) j++; 586674ae819SStefano Zampini total_counts=total_counts-j; 587674ae819SStefano Zampini if (j<temp_constraints) { 588674ae819SStefano Zampini for (k=j;k<Bt;k++) { 589674ae819SStefano Zampini singular_vals[k]=1.0/PetscSqrtReal(singular_vals[k]); 590674ae819SStefano Zampini } 591674ae819SStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 592674ae819SStefano Zampini BLASgemm_("N","N",&Bs,&Bt,&Bt,&one,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,correlation_mat,&Bt,&zero,temp_basis,&Bs); 593674ae819SStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 594674ae819SStefano Zampini /* copy POD basis into used quadrature memory */ 595674ae819SStefano Zampini for (k=0;k<Bt-j;k++) { 596674ae819SStefano Zampini for (ii=0;ii<size_of_constraint;ii++) { 597674ae819SStefano Zampini temp_quadrature_constraint[temp_indices[temp_start_ptr+k]+ii]=singular_vals[Bt-1-k]*temp_basis[(Bt-1-k)*size_of_constraint+ii]; 598674ae819SStefano Zampini } 599674ae819SStefano Zampini } 600674ae819SStefano Zampini } 601674ae819SStefano Zampini 602674ae819SStefano Zampini #else /* on missing GESVD */ 603674ae819SStefano Zampini PetscInt min_n = temp_constraints; 604674ae819SStefano Zampini if (min_n > size_of_constraint) min_n = size_of_constraint; 605674ae819SStefano Zampini dummy_int = Bs; 606674ae819SStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 607674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX) 608674ae819SStefano Zampini LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,singular_vals, 609674ae819SStefano Zampini &dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,work,&lwork,&lierr); 610674ae819SStefano Zampini #else 611674ae819SStefano Zampini LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,singular_vals, 612674ae819SStefano Zampini &dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,work,&lwork,rwork,&lierr); 613674ae819SStefano Zampini #endif 614674ae819SStefano Zampini if (lierr) { 615674ae819SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SVD Lapack routine %d",(int)lierr); 616674ae819SStefano Zampini } 617674ae819SStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 618674ae819SStefano Zampini /* retain eigenvalues greater than tol: note that lapack SVD gives eigs in descending order */ 619674ae819SStefano Zampini j = 0; 620674ae819SStefano Zampini while (j < min_n && singular_vals[min_n-j-1] < tol) j++; 621674ae819SStefano Zampini total_counts = total_counts-(PetscInt)Bt+(min_n-j); 622674ae819SStefano Zampini #endif 623674ae819SStefano Zampini } 624674ae819SStefano Zampini } 625674ae819SStefano Zampini /* free index sets of faces, edges and vertices */ 626674ae819SStefano Zampini for (i=0;i<n_ISForFaces;i++) { 627674ae819SStefano Zampini ierr = ISDestroy(&ISForFaces[i]);CHKERRQ(ierr); 628674ae819SStefano Zampini } 629674ae819SStefano Zampini ierr = PetscFree(ISForFaces);CHKERRQ(ierr); 630674ae819SStefano Zampini for (i=0;i<n_ISForEdges;i++) { 631674ae819SStefano Zampini ierr = ISDestroy(&ISForEdges[i]);CHKERRQ(ierr); 632674ae819SStefano Zampini } 633674ae819SStefano Zampini ierr = PetscFree(ISForEdges);CHKERRQ(ierr); 634674ae819SStefano Zampini ierr = ISDestroy(&ISForVertices);CHKERRQ(ierr); 635674ae819SStefano Zampini 636674ae819SStefano Zampini /* set quantities in pcbddc data structure */ 637674ae819SStefano Zampini /* n_vertices defines the number of point primal dofs */ 638674ae819SStefano Zampini /* n_constraints defines the number of averages (they can be point primal dofs if change of basis is requested) */ 639674ae819SStefano Zampini local_primal_size = total_counts; 640674ae819SStefano Zampini pcbddc->n_vertices = n_vertices; 641674ae819SStefano Zampini pcbddc->n_constraints = total_counts-n_vertices; 642674ae819SStefano Zampini pcbddc->local_primal_size = local_primal_size; 643674ae819SStefano Zampini 644674ae819SStefano Zampini /* Create constraint matrix */ 645674ae819SStefano Zampini /* The constraint matrix is used to compute the l2g map of primal dofs */ 646674ae819SStefano Zampini /* so we need to set it up properly either with or without change of basis */ 647674ae819SStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ConstraintMatrix);CHKERRQ(ierr); 648674ae819SStefano Zampini ierr = MatSetType(pcbddc->ConstraintMatrix,impMatType);CHKERRQ(ierr); 649674ae819SStefano Zampini ierr = MatSetSizes(pcbddc->ConstraintMatrix,local_primal_size,pcis->n,local_primal_size,pcis->n);CHKERRQ(ierr); 650674ae819SStefano Zampini /* compute a local numbering of constraints : vertices first then constraints */ 651674ae819SStefano Zampini ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr); 652674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr); 653674ae819SStefano Zampini ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_numbering);CHKERRQ(ierr); 654674ae819SStefano Zampini ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_permutation);CHKERRQ(ierr); 655674ae819SStefano Zampini total_counts=0; 656674ae819SStefano Zampini /* find vertices: subdomain corners plus dofs with basis changed */ 657674ae819SStefano Zampini for (i=0;i<local_primal_size;i++) { 658674ae819SStefano Zampini size_of_constraint=temp_indices[i+1]-temp_indices[i]; 659674ae819SStefano Zampini if (change_basis[i] || size_of_constraint == 1) { 660674ae819SStefano Zampini k=0; 661674ae819SStefano Zampini while(k < size_of_constraint && array_vector[temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1]] != 0.0) { 662674ae819SStefano Zampini k=k+1; 663674ae819SStefano Zampini } 664674ae819SStefano Zampini j=temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1]; 665674ae819SStefano Zampini array_vector[j] = 1.0; 666674ae819SStefano Zampini aux_primal_numbering[total_counts]=j; 667674ae819SStefano Zampini aux_primal_permutation[total_counts]=total_counts; 668674ae819SStefano Zampini total_counts++; 669674ae819SStefano Zampini } 670674ae819SStefano Zampini } 671674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr); 672674ae819SStefano Zampini /* permute indices in order to have a sorted set of vertices */ 673674ae819SStefano Zampini ierr = PetscSortIntWithPermutation(total_counts,aux_primal_numbering,aux_primal_permutation); 674674ae819SStefano Zampini /* nonzero structure */ 675674ae819SStefano Zampini ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 676674ae819SStefano Zampini for (i=0;i<total_counts;i++) { 677674ae819SStefano Zampini nnz[i]=1; 678674ae819SStefano Zampini } 679674ae819SStefano Zampini j=total_counts; 680674ae819SStefano Zampini for (i=n_vertices;i<local_primal_size;i++) { 681674ae819SStefano Zampini if (!change_basis[i]) { 682674ae819SStefano Zampini nnz[j]=temp_indices[i+1]-temp_indices[i]; 683674ae819SStefano Zampini j++; 684674ae819SStefano Zampini } 685674ae819SStefano Zampini } 686674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(pcbddc->ConstraintMatrix,0,nnz);CHKERRQ(ierr); 687674ae819SStefano Zampini ierr = PetscFree(nnz);CHKERRQ(ierr); 688674ae819SStefano Zampini /* set values in constraint matrix */ 689674ae819SStefano Zampini for (i=0;i<total_counts;i++) { 690674ae819SStefano Zampini j = aux_primal_permutation[i]; 691674ae819SStefano Zampini k = aux_primal_numbering[j]; 692674ae819SStefano Zampini ierr = MatSetValue(pcbddc->ConstraintMatrix,i,k,1.0,INSERT_VALUES);CHKERRQ(ierr); 693674ae819SStefano Zampini } 694674ae819SStefano Zampini for (i=n_vertices;i<local_primal_size;i++) { 695674ae819SStefano Zampini if (!change_basis[i]) { 696674ae819SStefano Zampini size_of_constraint=temp_indices[i+1]-temp_indices[i]; 697674ae819SStefano Zampini 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); 698674ae819SStefano Zampini total_counts++; 699674ae819SStefano Zampini } 700674ae819SStefano Zampini } 701674ae819SStefano Zampini ierr = PetscFree(aux_primal_numbering);CHKERRQ(ierr); 702674ae819SStefano Zampini ierr = PetscFree(aux_primal_permutation);CHKERRQ(ierr); 703674ae819SStefano Zampini /* assembling */ 704674ae819SStefano Zampini ierr = MatAssemblyBegin(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 705674ae819SStefano Zampini ierr = MatAssemblyEnd(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 706674ae819SStefano Zampini 707674ae819SStefano Zampini /* Create matrix for change of basis. We don't need it in case pcbddc->use_change_of_basis is FALSE */ 708674ae819SStefano Zampini if (pcbddc->use_change_of_basis) { 709674ae819SStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr); 710674ae819SStefano Zampini ierr = MatSetType(pcbddc->ChangeOfBasisMatrix,impMatType);CHKERRQ(ierr); 711674ae819SStefano Zampini ierr = MatSetSizes(pcbddc->ChangeOfBasisMatrix,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);CHKERRQ(ierr); 712674ae819SStefano Zampini /* work arrays */ 713674ae819SStefano Zampini /* we need to reuse these arrays, so we free them */ 714674ae819SStefano Zampini ierr = PetscFree(temp_basis);CHKERRQ(ierr); 715674ae819SStefano Zampini ierr = PetscFree(work);CHKERRQ(ierr); 716674ae819SStefano Zampini ierr = PetscMalloc(pcis->n_B*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 717674ae819SStefano Zampini ierr = PetscMalloc((nnsp_addone+nnsp_size)*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr); 718674ae819SStefano Zampini ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscScalar),&work);CHKERRQ(ierr); 719674ae819SStefano Zampini ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscBLASInt),&ipiv);CHKERRQ(ierr); 720674ae819SStefano Zampini for (i=0;i<pcis->n_B;i++) { 721674ae819SStefano Zampini nnz[i]=1; 722674ae819SStefano Zampini } 723674ae819SStefano Zampini /* Overestimated nonzeros per row */ 724674ae819SStefano Zampini k=1; 725674ae819SStefano Zampini for (i=pcbddc->n_vertices;i<local_primal_size;i++) { 726674ae819SStefano Zampini if (change_basis[i]) { 727674ae819SStefano Zampini size_of_constraint = temp_indices[i+1]-temp_indices[i]; 728674ae819SStefano Zampini if (k < size_of_constraint) { 729674ae819SStefano Zampini k = size_of_constraint; 730674ae819SStefano Zampini } 731674ae819SStefano Zampini for (j=0;j<size_of_constraint;j++) { 732674ae819SStefano Zampini nnz[temp_indices_to_constraint_B[temp_indices[i]+j]] = size_of_constraint; 733674ae819SStefano Zampini } 734674ae819SStefano Zampini } 735674ae819SStefano Zampini } 736674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(pcbddc->ChangeOfBasisMatrix,0,nnz);CHKERRQ(ierr); 737674ae819SStefano Zampini ierr = PetscFree(nnz);CHKERRQ(ierr); 738674ae819SStefano Zampini /* Temporary array to store indices */ 739674ae819SStefano Zampini ierr = PetscMalloc(k*sizeof(PetscInt),&is_indices);CHKERRQ(ierr); 740674ae819SStefano Zampini /* Set initial identity in the matrix */ 741674ae819SStefano Zampini for (i=0;i<pcis->n_B;i++) { 742674ae819SStefano Zampini ierr = MatSetValue(pcbddc->ChangeOfBasisMatrix,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr); 743674ae819SStefano Zampini } 744674ae819SStefano Zampini /* Now we loop on the constraints which need a change of basis */ 745674ae819SStefano Zampini /* Change of basis matrix is evaluated as the FIRST APPROACH in */ 746674ae819SStefano Zampini /* Klawonn and Widlund, Dual-primal FETI-DP methods for linear elasticity, (6.2.1) */ 747674ae819SStefano Zampini temp_constraints = 0; 748674ae819SStefano Zampini if (pcbddc->n_vertices < local_primal_size) { 749674ae819SStefano Zampini temp_start_ptr = temp_indices_to_constraint_B[temp_indices[pcbddc->n_vertices]]; 750674ae819SStefano Zampini } 751674ae819SStefano Zampini for (i=pcbddc->n_vertices;i<local_primal_size;i++) { 752674ae819SStefano Zampini if (change_basis[i]) { 753674ae819SStefano Zampini compute_submatrix = PETSC_FALSE; 754674ae819SStefano Zampini useksp = PETSC_FALSE; 755674ae819SStefano Zampini if (temp_start_ptr == temp_indices_to_constraint_B[temp_indices[i]]) { 756674ae819SStefano Zampini temp_constraints++; 757674ae819SStefano Zampini if (i == local_primal_size -1 || temp_start_ptr != temp_indices_to_constraint_B[temp_indices[i+1]]) { 758674ae819SStefano Zampini compute_submatrix = PETSC_TRUE; 759674ae819SStefano Zampini } 760674ae819SStefano Zampini } 761674ae819SStefano Zampini if (compute_submatrix) { 762674ae819SStefano Zampini if (temp_constraints > 1 || pcbddc->use_nnsp_true) { 763674ae819SStefano Zampini useksp = PETSC_TRUE; 764674ae819SStefano Zampini } 765674ae819SStefano Zampini size_of_constraint = temp_indices[i+1]-temp_indices[i]; 766674ae819SStefano Zampini if (useksp) { /* experimental TODO: reuse KSP and MAT instead of creating them each time */ 767674ae819SStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&temp_mat);CHKERRQ(ierr); 768674ae819SStefano Zampini ierr = MatSetType(temp_mat,impMatType);CHKERRQ(ierr); 769674ae819SStefano Zampini ierr = MatSetSizes(temp_mat,size_of_constraint,size_of_constraint,size_of_constraint,size_of_constraint);CHKERRQ(ierr); 770674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(temp_mat,size_of_constraint,NULL);CHKERRQ(ierr); 771674ae819SStefano Zampini } 772674ae819SStefano Zampini /* First _size_of_constraint-temp_constraints_ columns */ 773674ae819SStefano Zampini dual_dofs = size_of_constraint-temp_constraints; 774674ae819SStefano Zampini start_constraint = i+1-temp_constraints; 775674ae819SStefano Zampini for (s=0;s<dual_dofs;s++) { 776674ae819SStefano Zampini is_indices[0] = s; 777674ae819SStefano Zampini for (j=0;j<temp_constraints;j++) { 778674ae819SStefano Zampini for (k=0;k<temp_constraints;k++) { 779674ae819SStefano Zampini temp_basis[j*temp_constraints+k]=temp_quadrature_constraint[temp_indices[start_constraint+k]+s+j+1]; 780674ae819SStefano Zampini } 781674ae819SStefano Zampini work[j]=-temp_quadrature_constraint[temp_indices[start_constraint+j]+s]; 782674ae819SStefano Zampini is_indices[j+1]=s+j+1; 783674ae819SStefano Zampini } 784674ae819SStefano Zampini Bt = temp_constraints; 785674ae819SStefano Zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 786674ae819SStefano Zampini LAPACKgesv_(&Bt,&Bone,temp_basis,&Bt,ipiv,work,&Bt,&lierr); 787674ae819SStefano Zampini if ( lierr ) { 788674ae819SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GESV Lapack routine %d",(int)lierr); 789674ae819SStefano Zampini } 790674ae819SStefano Zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 791674ae819SStefano Zampini j = temp_indices_to_constraint_B[temp_indices[start_constraint]+s]; 792674ae819SStefano Zampini ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,temp_constraints,&temp_indices_to_constraint_B[temp_indices[start_constraint]+s+1],1,&j,work,INSERT_VALUES);CHKERRQ(ierr); 793674ae819SStefano Zampini if (useksp) { 794674ae819SStefano Zampini /* temp mat with transposed rows and columns */ 795674ae819SStefano Zampini ierr = MatSetValues(temp_mat,1,&s,temp_constraints,&is_indices[1],work,INSERT_VALUES);CHKERRQ(ierr); 796674ae819SStefano Zampini ierr = MatSetValue(temp_mat,is_indices[0],is_indices[0],1.0,INSERT_VALUES);CHKERRQ(ierr); 797674ae819SStefano Zampini } 798674ae819SStefano Zampini } 799674ae819SStefano Zampini if (useksp) { 800674ae819SStefano Zampini /* last rows of temp_mat */ 801674ae819SStefano Zampini for (j=0;j<size_of_constraint;j++) { 802674ae819SStefano Zampini is_indices[j] = j; 803674ae819SStefano Zampini } 804674ae819SStefano Zampini for (s=0;s<temp_constraints;s++) { 805674ae819SStefano Zampini k = s + dual_dofs; 806674ae819SStefano Zampini ierr = MatSetValues(temp_mat,1,&k,size_of_constraint,is_indices,&temp_quadrature_constraint[temp_indices[start_constraint+s]],INSERT_VALUES);CHKERRQ(ierr); 807674ae819SStefano Zampini } 808674ae819SStefano Zampini ierr = MatAssemblyBegin(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 809674ae819SStefano Zampini ierr = MatAssemblyEnd(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 810674ae819SStefano Zampini ierr = MatGetVecs(temp_mat,&temp_vec,NULL);CHKERRQ(ierr); 811674ae819SStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&temp_ksp);CHKERRQ(ierr); 812674ae819SStefano Zampini ierr = KSPSetOperators(temp_ksp,temp_mat,temp_mat,SAME_PRECONDITIONER);CHKERRQ(ierr); 813674ae819SStefano Zampini ierr = KSPSetType(temp_ksp,KSPPREONLY);CHKERRQ(ierr); 814674ae819SStefano Zampini ierr = KSPGetPC(temp_ksp,&temp_pc);CHKERRQ(ierr); 815674ae819SStefano Zampini ierr = PCSetType(temp_pc,PCLU);CHKERRQ(ierr); 816674ae819SStefano Zampini ierr = KSPSetUp(temp_ksp);CHKERRQ(ierr); 817674ae819SStefano Zampini for (s=0;s<temp_constraints;s++) { 818674ae819SStefano Zampini ierr = VecSet(temp_vec,0.0);CHKERRQ(ierr); 819674ae819SStefano Zampini ierr = VecSetValue(temp_vec,s+dual_dofs,1.0,INSERT_VALUES);CHKERRQ(ierr); 820674ae819SStefano Zampini ierr = VecAssemblyBegin(temp_vec);CHKERRQ(ierr); 821674ae819SStefano Zampini ierr = VecAssemblyEnd(temp_vec);CHKERRQ(ierr); 822674ae819SStefano Zampini ierr = KSPSolve(temp_ksp,temp_vec,temp_vec);CHKERRQ(ierr); 823674ae819SStefano Zampini ierr = VecGetArray(temp_vec,&array_vector);CHKERRQ(ierr); 824674ae819SStefano Zampini j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1]; 825674ae819SStefano Zampini /* last columns of change of basis matrix associated to new primal dofs */ 826674ae819SStefano Zampini ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,&temp_indices_to_constraint_B[temp_indices[start_constraint+s]],1,&j,array_vector,INSERT_VALUES);CHKERRQ(ierr); 827674ae819SStefano Zampini ierr = VecRestoreArray(temp_vec,&array_vector);CHKERRQ(ierr); 828674ae819SStefano Zampini } 829674ae819SStefano Zampini ierr = MatDestroy(&temp_mat);CHKERRQ(ierr); 830674ae819SStefano Zampini ierr = KSPDestroy(&temp_ksp);CHKERRQ(ierr); 831674ae819SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 832674ae819SStefano Zampini } else { 833674ae819SStefano Zampini /* last columns of change of basis matrix associated to new primal dofs */ 834674ae819SStefano Zampini for (s=0;s<temp_constraints;s++) { 835674ae819SStefano Zampini j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1]; 836674ae819SStefano Zampini ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,&temp_indices_to_constraint_B[temp_indices[start_constraint+s]],1,&j,&temp_quadrature_constraint[temp_indices[start_constraint+s]],INSERT_VALUES);CHKERRQ(ierr); 837674ae819SStefano Zampini } 838674ae819SStefano Zampini } 839674ae819SStefano Zampini /* prepare for the next cycle */ 840674ae819SStefano Zampini temp_constraints = 0; 841674ae819SStefano Zampini if (i != local_primal_size -1 ) { 842674ae819SStefano Zampini temp_start_ptr = temp_indices_to_constraint_B[temp_indices[i+1]]; 843674ae819SStefano Zampini } 844674ae819SStefano Zampini } 845674ae819SStefano Zampini } 846674ae819SStefano Zampini } 847674ae819SStefano Zampini /* assembling */ 848674ae819SStefano Zampini ierr = MatAssemblyBegin(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 849674ae819SStefano Zampini ierr = MatAssemblyEnd(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 850674ae819SStefano Zampini ierr = PetscFree(ipiv);CHKERRQ(ierr); 851674ae819SStefano Zampini ierr = PetscFree(is_indices);CHKERRQ(ierr); 852674ae819SStefano Zampini } 853674ae819SStefano Zampini /* free workspace no longer needed */ 854674ae819SStefano Zampini ierr = PetscFree(rwork);CHKERRQ(ierr); 855674ae819SStefano Zampini ierr = PetscFree(work);CHKERRQ(ierr); 856674ae819SStefano Zampini ierr = PetscFree(temp_basis);CHKERRQ(ierr); 857674ae819SStefano Zampini ierr = PetscFree(singular_vals);CHKERRQ(ierr); 858674ae819SStefano Zampini ierr = PetscFree(correlation_mat);CHKERRQ(ierr); 859674ae819SStefano Zampini ierr = PetscFree(temp_indices);CHKERRQ(ierr); 860674ae819SStefano Zampini ierr = PetscFree(change_basis);CHKERRQ(ierr); 861674ae819SStefano Zampini ierr = PetscFree(temp_indices_to_constraint);CHKERRQ(ierr); 862674ae819SStefano Zampini ierr = PetscFree(temp_indices_to_constraint_B);CHKERRQ(ierr); 863674ae819SStefano Zampini ierr = PetscFree(local_to_B);CHKERRQ(ierr); 864674ae819SStefano Zampini ierr = PetscFree(temp_quadrature_constraint);CHKERRQ(ierr); 865674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD) 866674ae819SStefano Zampini ierr = PetscFree(iwork);CHKERRQ(ierr); 867674ae819SStefano Zampini ierr = PetscFree(ifail);CHKERRQ(ierr); 868674ae819SStefano Zampini ierr = PetscFree(singular_vectors);CHKERRQ(ierr); 869674ae819SStefano Zampini #endif 870674ae819SStefano Zampini for (k=0;k<nnsp_size;k++) { 871674ae819SStefano Zampini ierr = VecDestroy(&localnearnullsp[k]);CHKERRQ(ierr); 872674ae819SStefano Zampini } 873674ae819SStefano Zampini ierr = PetscFree(localnearnullsp);CHKERRQ(ierr); 874674ae819SStefano Zampini PetscFunctionReturn(0); 875674ae819SStefano Zampini } 876674ae819SStefano Zampini 877674ae819SStefano Zampini #undef __FUNCT__ 878674ae819SStefano Zampini #define __FUNCT__ "PCBDDCAnalyzeInterface" 879674ae819SStefano Zampini PetscErrorCode PCBDDCAnalyzeInterface(PC pc) 880674ae819SStefano Zampini { 881674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 882674ae819SStefano Zampini PC_IS *pcis = (PC_IS*)pc->data; 883674ae819SStefano Zampini Mat_IS *matis = (Mat_IS*)pc->pmat->data; 884674ae819SStefano Zampini PetscInt bs,ierr,i,vertex_size; 885674ae819SStefano Zampini PetscViewer viewer=pcbddc->dbg_viewer; 886674ae819SStefano Zampini 887674ae819SStefano Zampini PetscFunctionBegin; 888674ae819SStefano Zampini /* Init local Graph struct */ 889674ae819SStefano Zampini ierr = PCBDDCGraphInit(pcbddc->mat_graph,matis->mapping);CHKERRQ(ierr); 890674ae819SStefano Zampini 891674ae819SStefano Zampini /* Set default CSR adjacency of local dofs if not provided by the user with PCBDDCSetLocalAdjacencyGraph */ 892674ae819SStefano Zampini if (!pcbddc->mat_graph->xadj || !pcbddc->mat_graph->adjncy) { 893674ae819SStefano Zampini Mat mat_adj; 894674ae819SStefano Zampini const PetscInt *xadj,*adjncy; 895674ae819SStefano Zampini PetscBool flg_row=PETSC_TRUE; 896674ae819SStefano Zampini 897674ae819SStefano Zampini ierr = MatConvert(matis->A,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr); 898674ae819SStefano Zampini ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr); 899674ae819SStefano Zampini if (!flg_row) { 900674ae819SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__); 901674ae819SStefano Zampini } 902674ae819SStefano Zampini ierr = PCBDDCSetLocalAdjacencyGraph(pc,i,xadj,adjncy,PETSC_COPY_VALUES);CHKERRQ(ierr); 903674ae819SStefano Zampini ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr); 904674ae819SStefano Zampini if (!flg_row) { 905674ae819SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 906674ae819SStefano Zampini } 907674ae819SStefano Zampini ierr = MatDestroy(&mat_adj);CHKERRQ(ierr); 908674ae819SStefano Zampini } 909674ae819SStefano Zampini 910674ae819SStefano Zampini /* Set default dofs' splitting if no information has been provided by the user with PCBDDCSetDofsSplitting */ 911674ae819SStefano Zampini vertex_size = 1; 912674ae819SStefano Zampini if (!pcbddc->n_ISForDofs) { 913674ae819SStefano Zampini IS *custom_ISForDofs; 914674ae819SStefano Zampini 915674ae819SStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 916674ae819SStefano Zampini ierr = PetscMalloc(bs*sizeof(IS),&custom_ISForDofs);CHKERRQ(ierr); 917674ae819SStefano Zampini for (i=0;i<bs;i++) { 918674ae819SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,pcis->n/bs,i,bs,&custom_ISForDofs[i]);CHKERRQ(ierr); 919674ae819SStefano Zampini } 920674ae819SStefano Zampini ierr = PCBDDCSetDofsSplitting(pc,bs,custom_ISForDofs);CHKERRQ(ierr); 921674ae819SStefano Zampini /* remove my references to IS objects */ 922674ae819SStefano Zampini for (i=0;i<bs;i++) { 923674ae819SStefano Zampini ierr = ISDestroy(&custom_ISForDofs[i]);CHKERRQ(ierr); 924674ae819SStefano Zampini } 925674ae819SStefano Zampini ierr = PetscFree(custom_ISForDofs);CHKERRQ(ierr); 926674ae819SStefano Zampini } else { /* mat block size as vertex size (used for elasticity) */ 927674ae819SStefano Zampini ierr = MatGetBlockSize(matis->A,&vertex_size);CHKERRQ(ierr); 928674ae819SStefano Zampini } 929674ae819SStefano Zampini 930674ae819SStefano Zampini /* Setup of Graph */ 931674ae819SStefano Zampini ierr = PCBDDCGraphSetUp(pcbddc->mat_graph,vertex_size,pcbddc->NeumannBoundaries,pcbddc->DirichletBoundaries,pcbddc->n_ISForDofs,pcbddc->ISForDofs,pcbddc->user_primal_vertices); 932674ae819SStefano Zampini 933674ae819SStefano Zampini /* Graph's connected components analysis */ 934674ae819SStefano Zampini ierr = PCBDDCGraphComputeConnectedComponents(pcbddc->mat_graph);CHKERRQ(ierr); 935674ae819SStefano Zampini 936674ae819SStefano Zampini /* print some info to stdout */ 937674ae819SStefano Zampini if (pcbddc->dbg_flag) { 938e49050b4SStefano Zampini ierr = PCBDDCGraphASCIIView(pcbddc->mat_graph,pcbddc->dbg_flag,viewer); 939674ae819SStefano Zampini } 940674ae819SStefano Zampini PetscFunctionReturn(0); 941674ae819SStefano Zampini } 942674ae819SStefano Zampini 943674ae819SStefano Zampini #undef __FUNCT__ 944674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGetPrimalVerticesLocalIdx" 945674ae819SStefano Zampini PetscErrorCode PCBDDCGetPrimalVerticesLocalIdx(PC pc, PetscInt *n_vertices, PetscInt *vertices_idx[]) 946674ae819SStefano Zampini { 947674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 948674ae819SStefano Zampini PetscInt *vertices,*row_cmat_indices,n,i,size_of_constraint,local_primal_size; 949674ae819SStefano Zampini PetscErrorCode ierr; 950674ae819SStefano Zampini 951674ae819SStefano Zampini PetscFunctionBegin; 952674ae819SStefano Zampini n = 0; 953674ae819SStefano Zampini vertices = 0; 954674ae819SStefano Zampini if (pcbddc->ConstraintMatrix) { 955674ae819SStefano Zampini ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&i);CHKERRQ(ierr); 956*b120a5c6SStefano Zampini for (i=0;i<local_primal_size;i++) { 957*b120a5c6SStefano Zampini ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 958*b120a5c6SStefano Zampini if (size_of_constraint == 1) n++; 959*b120a5c6SStefano Zampini ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 960*b120a5c6SStefano Zampini } 961*b120a5c6SStefano Zampini ierr = PetscMalloc(n*sizeof(PetscInt),&vertices);CHKERRQ(ierr); 962*b120a5c6SStefano Zampini n = 0; 963674ae819SStefano Zampini for (i=0;i<local_primal_size;i++) { 964674ae819SStefano Zampini ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 965674ae819SStefano Zampini if (size_of_constraint == 1) { 966674ae819SStefano Zampini vertices[n++]=row_cmat_indices[0]; 967674ae819SStefano Zampini } 968674ae819SStefano Zampini ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 969674ae819SStefano Zampini } 970674ae819SStefano Zampini } 971674ae819SStefano Zampini *n_vertices = n; 972674ae819SStefano Zampini *vertices_idx = vertices; 973674ae819SStefano Zampini PetscFunctionReturn(0); 974674ae819SStefano Zampini } 975674ae819SStefano Zampini 976674ae819SStefano Zampini #undef __FUNCT__ 977674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGetPrimalConstraintsLocalIdx" 978674ae819SStefano Zampini PetscErrorCode PCBDDCGetPrimalConstraintsLocalIdx(PC pc, PetscInt *n_constraints, PetscInt *constraints_idx[]) 979674ae819SStefano Zampini { 980674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 981674ae819SStefano Zampini PetscInt *constraints_index,*row_cmat_indices,*row_cmat_global_indices; 982674ae819SStefano Zampini PetscInt n,i,j,size_of_constraint,local_primal_size,local_size,max_size_of_constraint,min_index,min_loc; 983674ae819SStefano Zampini PetscBool *touched; 984674ae819SStefano Zampini PetscErrorCode ierr; 985674ae819SStefano Zampini 986674ae819SStefano Zampini PetscFunctionBegin; 987674ae819SStefano Zampini n = 0; 988674ae819SStefano Zampini constraints_index = 0; 989674ae819SStefano Zampini if (pcbddc->ConstraintMatrix) { 990674ae819SStefano Zampini ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&local_size);CHKERRQ(ierr); 991674ae819SStefano Zampini max_size_of_constraint = 0; 992674ae819SStefano Zampini for (i=0;i<local_primal_size;i++) { 993674ae819SStefano Zampini ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 994674ae819SStefano Zampini if (size_of_constraint > 1) { 995674ae819SStefano Zampini n++; 996674ae819SStefano Zampini } 997674ae819SStefano Zampini max_size_of_constraint = PetscMax(size_of_constraint,max_size_of_constraint); 998674ae819SStefano Zampini ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr); 999674ae819SStefano Zampini } 1000674ae819SStefano Zampini ierr = PetscMalloc(n*sizeof(PetscInt),&constraints_index);CHKERRQ(ierr); 1001674ae819SStefano Zampini ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&row_cmat_global_indices);CHKERRQ(ierr); 1002674ae819SStefano Zampini ierr = PetscMalloc(local_size*sizeof(PetscBool),&touched);CHKERRQ(ierr); 1003674ae819SStefano Zampini ierr = PetscMemzero(touched,local_size*sizeof(PetscBool));CHKERRQ(ierr); 1004674ae819SStefano Zampini n = 0; 1005674ae819SStefano Zampini for (i=0;i<local_primal_size;i++) { 1006674ae819SStefano Zampini ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1007674ae819SStefano Zampini if (size_of_constraint > 1) { 1008674ae819SStefano Zampini ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,row_cmat_indices,row_cmat_global_indices);CHKERRQ(ierr); 1009674ae819SStefano Zampini min_index = row_cmat_global_indices[0]; 1010674ae819SStefano Zampini min_loc = 0; 1011674ae819SStefano Zampini for (j=1;j<size_of_constraint;j++) { 1012674ae819SStefano Zampini /* there can be more than one constraint on a single connected component */ 1013674ae819SStefano Zampini if (min_index > row_cmat_global_indices[j] && !touched[row_cmat_indices[j]]) { 1014674ae819SStefano Zampini min_index = row_cmat_global_indices[j]; 1015674ae819SStefano Zampini min_loc = j; 1016674ae819SStefano Zampini } 1017674ae819SStefano Zampini } 1018674ae819SStefano Zampini touched[row_cmat_indices[min_loc]] = PETSC_TRUE; 1019674ae819SStefano Zampini constraints_index[n++] = row_cmat_indices[min_loc]; 1020674ae819SStefano Zampini } 1021674ae819SStefano Zampini ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr); 1022674ae819SStefano Zampini } 1023674ae819SStefano Zampini } 1024674ae819SStefano Zampini ierr = PetscFree(touched);CHKERRQ(ierr); 1025674ae819SStefano Zampini ierr = PetscFree(row_cmat_global_indices);CHKERRQ(ierr); 1026674ae819SStefano Zampini *n_constraints = n; 1027674ae819SStefano Zampini *constraints_idx = constraints_index; 1028674ae819SStefano Zampini PetscFunctionReturn(0); 1029674ae819SStefano Zampini } 1030674ae819SStefano Zampini 1031674ae819SStefano Zampini /* the next two functions has been adapted from pcis.c */ 1032674ae819SStefano Zampini #undef __FUNCT__ 1033674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplySchur" 1034674ae819SStefano Zampini PetscErrorCode PCBDDCApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 1035674ae819SStefano Zampini { 1036674ae819SStefano Zampini PetscErrorCode ierr; 1037674ae819SStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 1038674ae819SStefano Zampini 1039674ae819SStefano Zampini PetscFunctionBegin; 1040674ae819SStefano Zampini if (!vec2_B) { vec2_B = v; } 1041674ae819SStefano Zampini ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 1042674ae819SStefano Zampini ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr); 1043674ae819SStefano Zampini ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 1044674ae819SStefano Zampini ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr); 1045674ae819SStefano Zampini ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 1046674ae819SStefano Zampini PetscFunctionReturn(0); 1047674ae819SStefano Zampini } 1048674ae819SStefano Zampini 1049674ae819SStefano Zampini #undef __FUNCT__ 1050674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplySchurTranspose" 1051674ae819SStefano Zampini PetscErrorCode PCBDDCApplySchurTranspose(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 1052674ae819SStefano Zampini { 1053674ae819SStefano Zampini PetscErrorCode ierr; 1054674ae819SStefano Zampini PC_IS *pcis = (PC_IS*)(pc->data); 1055674ae819SStefano Zampini 1056674ae819SStefano Zampini PetscFunctionBegin; 1057674ae819SStefano Zampini if (!vec2_B) { vec2_B = v; } 1058674ae819SStefano Zampini ierr = MatMultTranspose(pcis->A_BB,v,vec1_B);CHKERRQ(ierr); 1059674ae819SStefano Zampini ierr = MatMultTranspose(pcis->A_BI,v,vec1_D);CHKERRQ(ierr); 1060674ae819SStefano Zampini ierr = KSPSolveTranspose(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr); 1061674ae819SStefano Zampini ierr = MatMultTranspose(pcis->A_IB,vec2_D,vec2_B);CHKERRQ(ierr); 1062674ae819SStefano Zampini ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr); 1063674ae819SStefano Zampini PetscFunctionReturn(0); 1064674ae819SStefano Zampini } 1065674ae819SStefano Zampini 1066674ae819SStefano Zampini #undef __FUNCT__ 1067674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSubsetNumbering" 1068674ae819SStefano Zampini 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[]) 1069674ae819SStefano Zampini { 1070674ae819SStefano Zampini Vec local_vec,global_vec; 1071674ae819SStefano Zampini IS seqis,paris; 1072674ae819SStefano Zampini VecScatter scatter_ctx; 1073674ae819SStefano Zampini PetscScalar *array; 1074674ae819SStefano Zampini PetscInt *temp_global_dofs; 1075674ae819SStefano Zampini PetscScalar globalsum; 1076674ae819SStefano Zampini PetscInt i,j,s; 1077674ae819SStefano Zampini PetscInt nlocals,first_index,old_index,max_local; 1078674ae819SStefano Zampini PetscMPIInt rank_prec_comm,size_prec_comm,max_global; 1079674ae819SStefano Zampini PetscMPIInt *dof_sizes,*dof_displs; 1080674ae819SStefano Zampini PetscBool first_found; 1081674ae819SStefano Zampini PetscErrorCode ierr; 1082674ae819SStefano Zampini 1083674ae819SStefano Zampini PetscFunctionBegin; 1084674ae819SStefano Zampini /* mpi buffers */ 1085674ae819SStefano Zampini MPI_Comm_size(comm,&size_prec_comm); 1086674ae819SStefano Zampini MPI_Comm_rank(comm,&rank_prec_comm); 1087674ae819SStefano Zampini j = ( !rank_prec_comm ? size_prec_comm : 0); 1088674ae819SStefano Zampini ierr = PetscMalloc(j*sizeof(*dof_sizes),&dof_sizes);CHKERRQ(ierr); 1089674ae819SStefano Zampini ierr = PetscMalloc(j*sizeof(*dof_displs),&dof_displs);CHKERRQ(ierr); 1090674ae819SStefano Zampini /* get maximum size of subset */ 1091674ae819SStefano Zampini ierr = PetscMalloc(n_local_dofs*sizeof(PetscInt),&temp_global_dofs);CHKERRQ(ierr); 1092674ae819SStefano Zampini ierr = ISLocalToGlobalMappingApply(l2gmap,n_local_dofs,local_dofs,temp_global_dofs);CHKERRQ(ierr); 1093674ae819SStefano Zampini max_local = 0; 1094674ae819SStefano Zampini if (n_local_dofs) { 1095674ae819SStefano Zampini max_local = temp_global_dofs[0]; 1096674ae819SStefano Zampini for (i=1;i<n_local_dofs;i++) { 1097674ae819SStefano Zampini if (max_local < temp_global_dofs[i] ) { 1098674ae819SStefano Zampini max_local = temp_global_dofs[i]; 1099674ae819SStefano Zampini } 1100674ae819SStefano Zampini } 1101674ae819SStefano Zampini } 1102674ae819SStefano Zampini ierr = MPI_Allreduce(&max_local,&max_global,1,MPIU_INT,MPI_MAX,comm); 1103674ae819SStefano Zampini max_global++; 1104674ae819SStefano Zampini max_local = 0; 1105674ae819SStefano Zampini if (n_local_dofs) { 1106674ae819SStefano Zampini max_local = local_dofs[0]; 1107674ae819SStefano Zampini for (i=1;i<n_local_dofs;i++) { 1108674ae819SStefano Zampini if (max_local < local_dofs[i] ) { 1109674ae819SStefano Zampini max_local = local_dofs[i]; 1110674ae819SStefano Zampini } 1111674ae819SStefano Zampini } 1112674ae819SStefano Zampini } 1113674ae819SStefano Zampini max_local++; 1114674ae819SStefano Zampini /* allocate workspace */ 1115674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr); 1116674ae819SStefano Zampini ierr = VecSetSizes(local_vec,PETSC_DECIDE,max_local);CHKERRQ(ierr); 1117674ae819SStefano Zampini ierr = VecSetType(local_vec,VECSEQ);CHKERRQ(ierr); 1118674ae819SStefano Zampini ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr); 1119674ae819SStefano Zampini ierr = VecSetSizes(global_vec,PETSC_DECIDE,max_global);CHKERRQ(ierr); 1120674ae819SStefano Zampini ierr = VecSetType(global_vec,VECMPI);CHKERRQ(ierr); 1121674ae819SStefano Zampini /* create scatter */ 1122674ae819SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,n_local_dofs,local_dofs,PETSC_COPY_VALUES,&seqis);CHKERRQ(ierr); 1123674ae819SStefano Zampini ierr = ISCreateGeneral(comm,n_local_dofs,temp_global_dofs,PETSC_COPY_VALUES,&paris);CHKERRQ(ierr); 1124674ae819SStefano Zampini ierr = VecScatterCreate(local_vec,seqis,global_vec,paris,&scatter_ctx);CHKERRQ(ierr); 1125674ae819SStefano Zampini ierr = ISDestroy(&seqis);CHKERRQ(ierr); 1126674ae819SStefano Zampini ierr = ISDestroy(&paris);CHKERRQ(ierr); 1127674ae819SStefano Zampini /* init array */ 1128674ae819SStefano Zampini ierr = VecSet(global_vec,0.0);CHKERRQ(ierr); 1129674ae819SStefano Zampini ierr = VecSet(local_vec,0.0);CHKERRQ(ierr); 1130674ae819SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 1131674ae819SStefano Zampini if (local_dofs_mult) { 1132674ae819SStefano Zampini for (i=0;i<n_local_dofs;i++) { 1133674ae819SStefano Zampini array[local_dofs[i]]=(PetscScalar)local_dofs_mult[i]; 1134674ae819SStefano Zampini } 1135674ae819SStefano Zampini } else { 1136674ae819SStefano Zampini for (i=0;i<n_local_dofs;i++) { 1137674ae819SStefano Zampini array[local_dofs[i]]=1.0; 1138674ae819SStefano Zampini } 1139674ae819SStefano Zampini } 1140674ae819SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 1141674ae819SStefano Zampini /* scatter into global vec and get total number of global dofs */ 1142674ae819SStefano Zampini ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1143674ae819SStefano Zampini ierr = VecScatterEnd (scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1144674ae819SStefano Zampini ierr = VecSum(global_vec,&globalsum);CHKERRQ(ierr); 1145674ae819SStefano Zampini *n_global_subset = (PetscInt)globalsum; 1146674ae819SStefano Zampini /* Fill global_vec with cumulative function for global numbering */ 1147674ae819SStefano Zampini ierr = VecGetArray(global_vec,&array);CHKERRQ(ierr); 1148674ae819SStefano Zampini ierr = VecGetLocalSize(global_vec,&s);CHKERRQ(ierr); 1149674ae819SStefano Zampini nlocals = 0; 1150674ae819SStefano Zampini first_index = -1; 1151674ae819SStefano Zampini first_found = PETSC_FALSE; 1152674ae819SStefano Zampini for (i=0;i<s;i++) { 1153674ae819SStefano Zampini if (!first_found && array[i] > 0.0) { 1154674ae819SStefano Zampini first_found = PETSC_TRUE; 1155674ae819SStefano Zampini first_index = i; 1156674ae819SStefano Zampini } 1157674ae819SStefano Zampini nlocals += (PetscInt)array[i]; 1158674ae819SStefano Zampini } 1159674ae819SStefano Zampini ierr = MPI_Gather(&nlocals,1,MPIU_INT,dof_sizes,1,MPIU_INT,0,comm);CHKERRQ(ierr); 1160674ae819SStefano Zampini if (!rank_prec_comm) { 1161674ae819SStefano Zampini dof_displs[0]=0; 1162674ae819SStefano Zampini for (i=1;i<size_prec_comm;i++) { 1163674ae819SStefano Zampini dof_displs[i] = dof_displs[i-1]+dof_sizes[i-1]; 1164674ae819SStefano Zampini } 1165674ae819SStefano Zampini } 1166674ae819SStefano Zampini ierr = MPI_Scatter(dof_displs,1,MPIU_INT,&nlocals,1,MPIU_INT,0,comm);CHKERRQ(ierr); 1167674ae819SStefano Zampini if (first_found) { 1168674ae819SStefano Zampini array[first_index] += (PetscScalar)nlocals; 1169674ae819SStefano Zampini old_index = first_index; 1170674ae819SStefano Zampini for (i=first_index+1;i<s;i++) { 1171674ae819SStefano Zampini if (array[i] > 0.0) { 1172674ae819SStefano Zampini array[i] += array[old_index]; 1173674ae819SStefano Zampini old_index = i; 1174674ae819SStefano Zampini } 1175674ae819SStefano Zampini } 1176674ae819SStefano Zampini } 1177674ae819SStefano Zampini ierr = VecRestoreArray(global_vec,&array);CHKERRQ(ierr); 1178674ae819SStefano Zampini ierr = VecSet(local_vec,0.0);CHKERRQ(ierr); 1179674ae819SStefano Zampini ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1180674ae819SStefano Zampini ierr = VecScatterEnd (scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1181674ae819SStefano Zampini /* get global ordering of local dofs */ 1182674ae819SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 1183674ae819SStefano Zampini if (local_dofs_mult) { 1184674ae819SStefano Zampini for (i=0;i<n_local_dofs;i++) { 1185674ae819SStefano Zampini temp_global_dofs[i] = (PetscInt)array[local_dofs[i]]-local_dofs_mult[i]; 1186674ae819SStefano Zampini } 1187674ae819SStefano Zampini } else { 1188674ae819SStefano Zampini for (i=0;i<n_local_dofs;i++) { 1189674ae819SStefano Zampini temp_global_dofs[i] = (PetscInt)array[local_dofs[i]]-1; 1190674ae819SStefano Zampini } 1191674ae819SStefano Zampini } 1192674ae819SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 1193674ae819SStefano Zampini /* free workspace */ 1194674ae819SStefano Zampini ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr); 1195674ae819SStefano Zampini ierr = VecDestroy(&local_vec);CHKERRQ(ierr); 1196674ae819SStefano Zampini ierr = VecDestroy(&global_vec);CHKERRQ(ierr); 1197674ae819SStefano Zampini ierr = PetscFree(dof_sizes);CHKERRQ(ierr); 1198674ae819SStefano Zampini ierr = PetscFree(dof_displs);CHKERRQ(ierr); 1199674ae819SStefano Zampini /* return pointer to global ordering of local dofs */ 1200674ae819SStefano Zampini *global_numbering_subset = temp_global_dofs; 1201674ae819SStefano Zampini PetscFunctionReturn(0); 1202674ae819SStefano Zampini } 1203